{ "cells": [ { "cell_type": "markdown", "source": [ "# Monitoring self-consistent field calculations\n", "\n", "The `self_consistent_field` function takes as the `callback`\n", "keyword argument one function to be called after each iteration.\n", "This function gets passed the complete internal state of the SCF\n", "solver and can thus be used both to monitor and debug the iterations\n", "as well as to quickly patch it with additional functionality.\n", "\n", "This example discusses a few aspects of the `callback` function\n", "taking again our favourite silicon example." ], "metadata": {} }, { "cell_type": "markdown", "source": [ "We setup silicon in an LDA model using the ASE interface\n", "to build the silicon lattice,\n", "see Creating slabs with ASE for more details." ], "metadata": {} }, { "outputs": [], "cell_type": "code", "source": [ "using DFTK\n", "using PyCall\n", "\n", "silicon = pyimport(\"ase.build\").bulk(\"Si\")\n", "atoms = [ElementPsp(el.symbol, psp=load_psp(\"hgh/lda/si-q4.hgh\"))\n", " for el in load_atoms(silicon)]\n", "positions = load_positions(silicon)\n", "lattice = load_lattice(silicon);\n", "\n", "model = model_LDA(lattice, atoms, positions)\n", "basis = PlaneWaveBasis(model; Ecut=5, kgrid=[3, 3, 3]);" ], "metadata": {}, "execution_count": 1 }, { "cell_type": "markdown", "source": [ "DFTK already defines a few callback functions for standard\n", "tasks. One example is the usual convergence table,\n", "which is defined in the callback `ScfDefaultCallback`.\n", "Another example is `ScfPlotTrace`, which records the total\n", "energy at each iteration and uses it to plot the convergence\n", "of the SCF graphically once it is converged.\n", "For details and other callbacks\n", "see [`src/scf/scf_callbacks.jl`](https://dftk.org/blob/master/src/scf/scf_callbacks.jl).\n", "\n", "!!! note \"Callbacks are not exported\"\n", " Callbacks are not exported from the DFTK namespace as of now,\n", " so you will need to use them, e.g., as `DFTK.ScfDefaultCallback`\n", " and `DFTK.ScfPlotTrace`." ], "metadata": {} }, { "cell_type": "markdown", "source": [ "In this example we define a custom callback, which plots\n", "the change in density at each SCF iteration after the SCF\n", "has finished. For this we first define the empty plot canvas\n", "and an empty container for all the density differences:" ], "metadata": {} }, { "outputs": [], "cell_type": "code", "source": [ "using Plots\n", "p = plot(yaxis=:log)\n", "density_differences = Float64[];" ], "metadata": {}, "execution_count": 2 }, { "cell_type": "markdown", "source": [ "The callback function itself gets passed a named tuple\n", "similar to the one returned by `self_consistent_field`,\n", "which contains the input and output density of the SCF step\n", "as `ρin` and `ρout`. Since the callback gets called\n", "both during the SCF iterations as well as after convergence\n", "just before `self_consistent_field` finishes we can both\n", "collect the data and initiate the plotting in one function." ], "metadata": {} }, { "outputs": [], "cell_type": "code", "source": [ "using LinearAlgebra\n", "\n", "function plot_callback(info)\n", " if info.stage == :finalize\n", " plot!(p, density_differences, label=\"|ρout - ρin|\", markershape=:x)\n", " else\n", " push!(density_differences, norm(info.ρout - info.ρin))\n", " end\n", " info\n", "end\n", "callback = DFTK.ScfDefaultCallback() ∘ plot_callback;" ], "metadata": {}, "execution_count": 3 }, { "cell_type": "markdown", "source": [ "Notice that for constructing the `callback` function we chained the `plot_callback`\n", "(which does the plotting) with the `ScfDefaultCallback`, such that when using\n", "the `plot_callback` function with `self_consistent_field` we still get the usual\n", "convergence table printed. We run the SCF with this callback ..." ], "metadata": {} }, { "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "n Energy log10(ΔE) log10(Δρ) α Diag\n", "--- --------------- --------- --------- ---- ----\n", " 1 -7.846733854537 -0.70 0.80 4.0\n", " 2 -7.852286547435 -2.26 -1.54 0.80 1.0\n", " 3 -7.852621368918 -3.48 -2.53 0.80 3.2\n", " 4 -7.852621946492 -6.24 -3.36 0.80 2.4\n", " 5 -7.852621955384 -8.05 -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+gvaeTAAAgAElEQVR4nO3deWBU1d248XPvPVlJwg4JYQ1E2SMhQQJBwr5mI0CplrgitRv+UCv4vlXeuhT0bd+mUhVsrSJYKxCCIWxhqwTCYgKEXTbLpuyyZL/3zu+PUVQECckkd2bu8/krGZKZL5bmyblzzozicDgEAAB2pVo9AAAAViKEAABbI4QAAFsjhAAAWyOEAABbI4QAAFsjhAAAWyOEAABbI4QAAFuroxDm5OQ8/fTTWVlZdfNwAABUUV2EMCsrKyMjY+zYse+8884HH3xwqy87efLkmTNnqn63hmG4YjpUh2mavDifhfjHby3TNK0ewb5q44dPXYRwzpw5L7zwQu/evV988cU5c+bc6sv+/Oc/z58/v+p3W1JS4orpUB0VFRX8LLYQ//gt5HA4+O9vobKyMpf/IlIXITxy5EhkZKQQIjIy8siRI3XwiAAAVFFdhPD6MlZRFC4pAADcSk1DeOrUqXXr1t3w3F5paWlWVtYHH3xw7tw5IUSbNm0+//xzIcSxY8fatm1bw0cEAMCFahTCqKiojh07jhw5csOGDddvvHbtWu/evV9//fXly5d36dLlwIEDDz300MyZMw8fPvzKK688+OCDNR0ZAADXkTX55iVLlrRt2zY6Ovq7N86bNy84ODg3N1dV1aeffvoPf/jDe++9ZxjGiy++2L9//0mTJtXkEc+Viab+N95YrAtFiMAa/VUAADZVo3pERET88Mbs7OwxY8aoqiqESEtLGzVqlBDiwQcfvO1asKKi4umnn3766aednyYnJ7///vs3fM2DG33HtjHGtTaEENeuXRNClBjKTzb6TO2oDwjl2ce6U1ZWJqWUkt8+rOH8xw9LOHeNcnzIKqWlpRUVFZqmVfHrAwMDb/vFrv9BdurUqfDwcOfHLVu2vHTpUnFxcb169W77jb6+vv/7v//71FNP/cjXLBoqUnN1Hz91YgdVCKEFBP8sV5/USU3q4OeS4VFFPj4+hNBawcHBVo9gUw6HQ1XVoKAgqwexKU3T/Pz8qh7CqnD9DzLDMK6P6PzAhQfOAqVYMkSm5upCiGFNlAdy9Yci1Qc68EJxAIBqcn1CwsLCzp496/z4zJkzQUFBISEhLrx/Zwv/8ZkZv9rnkbuoIACgRlxfkYSEhNWrVzs/Xr16dUJCgssfwiGE6RDnytRp242tZ7lSDwCovhpdGs3IyNi3b9/JkyfffvvtdevWPfXUU3fdddekSZNmz57961//Oiws7NVXX122bJmrZnUq0UVKrj7pbnVIk5Lh6/1GrtJHtVL/L05rzLOEAOrQxx9/fEevCokaGjFixMMPP1wb91yjEEZGRgYEBPTs2dP5qfPZ4+bNmxcUFMybN6+0tHTDhg333HOPC8b8Rokukr95XvDqVUdeok/iKv1MmaP7Yv0PsWp6JJdJAdSRTZs2+fn5JSUlWT2ILeTl5a1Zs8YdQzhy5Mib3h4eHj59+vSa3POt/GSdPuludXzE18ELlCJ7mEzO1f/rHvUve813PzPfitfuqq/UxkMDwA26des2btw4q6ewhcrKypycnFq6cw/b/v5uf3nDJdBAKZYOkYoQj3dU/7rPjM/Wf9FZfe4ezZfFIQCgCjwsFzd9IjBQigAppCqmdFW3Jstt5xwxWfrmM2yiAQDcnoeF8LbaBSvLh8mXY9SfrjfSNxgXyq0eCADg3rwthE6JrdWiMbKhn+i+WJ93iJdeAwDckneGUAhR31dkxGnZQ7XX95oDcvSDl7lSCsB673xmLjx2m9/OP/6P+eZ+foOvO14bQqfoJkp+kkxpo/bL1mcUGuUue603AKiOCRHq3APmB0du2bkVJxyziswJEV7+w9mteP9/a+cmmh2pcs8l0S1TX3eapSEAyzg3uv/js5u3cMUJx0s7jWVDZcNqvUJIWFjY1atXazqiEEIIXdfXrFlT8zfZuHbtWpMmTX78fjp37nz8+PEaPlBNeH8IncLrKYsGaX+8V334EyN9g3G+zOqBANiVs4XvHLyxhTWsoBCivLzcVe8PVVJSMmTIkJq/ZYKvr+/zzz+vKD92vNuFY1ePh50jrKHE1urAFuqLO4zOiypfitEmdVQ5ew+g7gVKsXSoTF6tCyHub68KV1TwBg6H44MPPti9e3dkZOTEiRN9fX2dt69atWrz5s2hoaH3339//fr1hRBr165t0aJFp06dhBDHjx/fuXNnUlJSVlaWEOJvf/ubqqopKSnNmjWryoNeu3Zt7ty5Fy9e7NOnj/MVVxRF8ff/+u3UMzMze/XqtWbNmqNHjw4ZMqRfv36u+avWmF1WhNfVk2JmrLZqhPzbQXNAjn7gK66UArBAPSmWDv16XejyCgohJk+e/O6773bo0GHRokVpaWnOG59//vmpU6eGhoYWFBT06NHj8uXLQoi333773//+t/MLioqKXn311eo9YklJSXx8/OHDh1u1avXkk0/+6U9/EkKUlpZOnjzZNE0hxMsvvzxhwoSjR48GBQUlJSVt2rTJBX9PV7DXivC6Ho2VzYnS+Uo0j3VU/yda83PluzwCsCPdFCNW6Xf0LaZDPJFnmELENFHGr7uz780aIuvd4kf44cOHMzMzjx8/HhgYOHHixDZt2mzbti0yMvK1117bvXt3hw4dhBDDhg2bM2fOb3/725veQ0pKihDiscceq/r7by9YsKBFixZvvPGGECI2NnbgwIG/+c1vbviaxMTEZ599Vghx9uzZpUuX9u3bt4p3XqtsGkLxzSaacRHKlHyzW6b+Rh9tcDgXSgFUn6qIZ7vf2e/U+WfMs2VCFaJfqHJf6J1dovO79Zfv3r27W7dugYGBQgg/P78ePXrs2bNHCNGwYUNnBYUQ/fr1c954R06fPu0smb+//9tvv/3dPyoqKurdu7fz4x49epSXl584caJx48bf/Zrrb8MQHh5eVFR0p49eS+wbQqcWgcrCQVr2cfOxjcZ9ocofe2tN/a2eCYBnUhVxR79PrzjhWHnKsXG09NVE8mq9YwOH8/nCmrt27VpZ2bd7AktLSwMDAwMDA0tLS6/fWFJS4iylpmnXN8WUlJT8+D2HhIQ4L7T+cKVYXFzsfA8iIYSu65WVlc77/y5Nc8eLb7Z7jvCmElur+8bKFvVEl0WVcw+YPG0IoLZ993nB7z5f6Kr7Lygo+Oyzz4QQx44dKygo6NOnT2RkZGBgYHZ2thCiuLh48eLFzjdOb9u2bWFhoRDC4XBkZmY6v71evXo+Pj4XLly44W6DgoJSUlJSUlJGjx79wwddunSps7UfffRRZGRk8+bNXfXXqVV2XxFeFyjFzFhtQoQ6Oc9YcNh8K17r1IArpQBqxQ93xzhbmJL77T7SGoqKinrggQfCwsK2bt36+9//vnXr1kKId9555+GHH37jjTcOHDiQkJAwfvx4IcSkSZP69u07YMCAy5cvd+7c2fntmqb94he/iIqKatmy5XvvvdelS5eqPGhERERcXFx4ePinn366cOHCmv8t6gYh/J57Giv5SfJvB81+2fpjHdUZ0Zq/O67jAXiwW+0RrSfF0iEy2UUtbN68+eLFi3fu3NmmTZvQ0FDnjcOGDTt8+PD+/ftDQ0PDw8OdN7Zt2/bgwYP79u1r06ZNo0aNrl9T/fOf//zqq68WFxcHBwdX8UETEhIefvjhw4cPd+7c2fldwcHBZ86cUVVVCLF+/fqAgADnVz7xxBM1P6ToKoTwRqoiHu+oJrZWn91mdFusv9FXG8ImGgAusvyEY+YuY/kwWd/3Jn8aKMWSwTJpte6viTFta9pCPz+/e++998aHCAzs2bPnDTcGBQX16tXL+bGPj8/12319fa8fQKyixo0bf3eDjKIo188ghoSEXL/9+uFCd0AIby4sUMxL0JYdd0zaaPRqqszuozULsHomAJ7v7voie+jNK+gU5CM+HipPl1Rnr0KrVq2ca6+uXbvW/Xpr5MiRVTx3f4Pw8PCqH9KoDYTwx4xurQxsIV8tMrplVj4Xpf26i8pL0QCoifYht/8hEuQj7qpfnZ81u3btcn7Qo0ePHj16VOMeamLs2LHV+8ZPPvnEtZPcKXaN3kagFDOitdUj5IdHzYQcfR+vRAMA3oUQVklUI2VTovxZB/W+bH3adqPMXZ7iBQDUFCGsKucmmt1pPqeLRdfF+upTLA0BwBsQwjvj3ESTEadNzjPGrzXOlt7+WwAA7ozNMtUxqpUyII1NNICtRUREvPjiix999JHVg9jCxYsX4+Pja+nOCWE1OTfRjGmr/jzPWHjMnBOvdWlIDAEbmTRpUkxMjNVT2EirVq1q6Z4JYY10b6RsSpLvHzIHLdd/EqG+Eqvd6l1RAHgZVVV/eDIdnojnCGtKESI9Ut05xudSuei+WF95kk00AOBJWL+4RmiAmJegrf/C8fM846764q2+Wng9rpQCgAdgRehKA8KUnamyZxMlOkvP2GPyfk4A4P4IoYsFSDEjWlszQn50zOy1VC84TwwBwK0RwlrRrZGSlyh/00UdvUqfkm9cq7R6IADALRDC2vLdTTRRmfqKEywNAcAdEcLa1TxAzEvQ/n6f9v+2GImr9ZPF5BAA3AshrAsJYcqOVNmzidIzS8/YYxrUEADcBiGsI85NNHmj5cfHzV5L9U/ZRAMA7oEQ1qnI+sqakXJKFzWRTTQA4B4IYV1zbqLZO9anzBCdFumLj5lWTwQAtkYIrdHIT8yJ1+YnaP9dYCau1k+wiQYALEIIrdQ/TCkaI+Obq/dk6rN2sYkGACxACC3mo4pno9QtSXLNaTM2S99+jhgCQJ0ihG4hsr6yeoR8squatFqfkm9cZRMNANQVQugunJto9n2ziWYRm2gAoE4QQvfS0E/Midc+GKA9X2AmrtaPX+NKKQDULkLoju4LVXaNkfHN1R5L2EQDALWLELop5yaarcly7WkzJkvfxiYaAKgdhNCtdQhRVo+Qz92jJq3WJ+cZV9hEAwCuRgg9wLh26v6xPv6a6LxIn3eITTQA4EqE0DM09BMZcdo/B2izdpmJq/X/sIkGAFyEEHqSfqHKzjEyvrkavUSfUWhUsjgEgBojhB7GuYlmW7LMP+OIzdK3nGVpCAA1Qgg9UvsQZdUI+V891NRcNtEAQI0QQg82rp26fxybaACgRgihZ2vgKzLitA8HaK8WmaNX6Z9f5UopANwZQugN4kOVHalySLjaa6k+o9CoYHEIAFVGCL2EjyqmdFW3JMutZx2xWXo+m2gAoGoIoVeJCFZWDJcvxag/WWukbzAulFs9EAC4PULohRJbq0VpsqGf6L6YTTQAcBuE0Ds5N9F8PFT7y15zYI7+2WWulALAzRFCb9azibIlSSa3UeOz2UQDADdHCL2cVMWUrurWZLntnCMmS998hqUhAHwPIbSFdsHK8mHy5Rj1p+vZRAMA30MIbSSxtVo0hk00APA9hNBe6vuKjDgte6j2+l5zQI5+kE00AGyPENpRdBMlP0mmtFH7ZeszCo1yw+qBAMA6hNCmnJtodqTKPZdEt0x93WmWhgBsihDaWng9ZdEg7Y/3qg9/YqRvMM6XWT0QANQ5QgiR2FrdN1a2qCc6L6qce8BkbQjAVgghhBCinhQzY7VVI+TfDpoDcvQDX1FDAHZBCPGtHo2VzYkytY0an61P284mGgC2QAjxPc5NNEVp8sgV0S1TX3OKpSEAL0cIcRMtApWFg7Q/3qs+ttFI32CcYxMNAO9FCHFL1zfRdFlUOfeAearE8eGRm7wezZv7zVK97qcDANeQVg8AtxYoxcxYbUKEOjnPmH9YSFW5VGE+0enb35/+61PjfJnw72ThjABQI4QQt3dPYyU/Sf7toPncduNEsVJpisfbC/FNBd+K1xSrJwSAaiOEqBJVEY93VBNbq89sM5771Dx4SQ3xdVysUKggAE9HCHEHwgLF/ARt4THz/g2OMD9x8CdUEIDHY7MM7tjOC44JbR2morT/l77nEucrAHg2Qog743xecG6c+VmaCA0UcR/rbx/grQ0BeDBCiDvw3d0xvqrYkiRjm6ozdpjpG4xrlVYPBwDVQghRVUevOirN7+0R9VXF8mHa0BZKsI+IydJ3XeQyKQDPw2YZVFVEsPJqL+2GG/018Y/+mhBi0TFz2Ap9epQ2pSu/XQHwJPzMgmuMbad+Mlq+e8hMW2NcKrd6GgCoMkIIl7mrvrIlSbasJ6Kz9C1nuUwKwDMQQriSnyYy4rT/660m5+qzdvEevwA8ACGE66W0Ubcny4+Pmym5xgUukwJwb4QQtaJ1kPLvUbJHYxG9RM/7kpUhAPdFCFFbpCpmRGt/76dNWG/MKDS4TgrAPRFC1K7B4cqWJG39aceQFfoXJVZPAwA/QAhR61rWU9aNkv1ClZ5ZlatOsjAE4F4IIeqCpogZ0dqCAfLRjcaUfKOSVycF4DYIIerOgDBlR6r87LJjyAr9VDFLQwBugRCiTjX1F8uHy9Q2auxSPecELQRgPUKIuqYIMaWrmjVE/mazMSXfqOAyKQBLEUJYo1dTpTBVni4RfT/Wj15laQjAMoQQlqnvKz4apP2sg9p7qf6voywMAViDEMJKzsukK4bL//7UTN9glOpWDwTAfgghrNeziVKQKisdom+2fugyl0kB1ClCCLcQ4iP+OUB7sqvab5m+4DCXSQHUHUIIN5Ieqa4bJWfuMtM3GMVcJgVQJwgh3EvnBsq2ZOkvRUyWXnSRy6QAah0hhNsJkGJuvDY9Sh20XM/Yw2VSALWLEMJNpUeqeYnyH5+ZY9caX1VYPQ0A70UI4b7urq9sSZbhgaLHEn3rWS6TAqgVhBBuzV8TGXHaH+9VU3L1jD28uS8A1yOE8ABj2qpbkuWHR83UXONiudXTAPAuhBCeoU2QsnG07NhA9FiibzrDyhCAyxBCeAypipmx2utxatoafUahwXVSAC5BCOFhktqo21Pk2tOOYSv1L0utngaA5yOE8Dyt6inrR8q+zZXoJZW5p1gYAqgRQgiPJFUxI1qbnyAf+cSYUWgY1BBAdRFCeLCBLZSCFLnlrGPwcv10CTEEUB2EEJ6tWYBYMVymtFFjsvQVJ2ghgDtGCOHxnO/u+88B8vE8Y0q+UcmrkwK4E4QQXqJ/mLIjVR654ojP1o9dZWkIoKoIIbxHE3+RPUze3169d6m+8BgLQwBVQgjhVZyXSXOGyenbzSn5Rrlh9UAA3B4hhBeKbaoUpsozpaJPtn74CpdJAfwYQgjvFOIjPhyoTemixmfr/zzCZVIAt0QI4c3SI9UVw+ULhWb6BqNEt3oaAG6JEMLL9WisFKRIU4iYLH3PJS6TArgRIYT3C/YR8xO0aVHqgBw9Yw+XSQF8DyGEXaRHqhtHy3c+M9M3GNcqrZ4GgNsghLCRjg2UrcmyoZ+IydJ3XuAyKQAhCCHsxl8TGXHayzHq8JVcJgUgBCGEPaW1U7ckyX8eNcesMS6VWz0NAEsRQthU22Dl36Nkq3qixxI9/yyXSQH7IoSwLz9NZMRpGXFqSq4+o9AwqSFgS4QQdpfcRt2eLNeccqTkGhe4TArYDyEEROsgZcMoGd1ERC/RN37JwhCwF0IICCGEVMWMaO3v/bSfrje4TArYCiEEvjU4XClIkZvPOAYv178osXoaAHWCEALf0zxArBgu7wtTopdUrjrJwhDwfoQQuJGmiBnR2j8Hykc3GlPyjUqO3QNejRACN5cQpuxIlYeuOPot0z+/ytIQ8FqEELilpv4iZ5j8aYTaJ1vPOUELAe9ECIEfowgxpauaNUT+ZrMxJd+o4DIp4HUIIXB7vZoqhanyixLR92P9yBWWhoBXIYRAldT3Ff8apP2sgxr3sf6voywMAe9BCIGqcl4mXTlc/venZvoGo0S3eiAArkAIgTsT3UQpTJW6Q/Raqu+9xGVSwOMRQuCOBfuIDwZov+2uDlquzz/MZVLAsxFCoJrSI9X1o+SsXWb6BqOYy6SAxyKEQPV1aqBsS5YN/URMll50kcukgEcihECNBEiREadNj1IHLdcz9nCZFPA8hBBwgfRINS9RvnvITFtjfFVh9TQA7gQhBFzj7vpKfpJsWU/0WKJvOctlUsBjEELAZfw1kRGn/eleNTVXz9jDm/sCnoEQAi6W2lbdkiz/ddRMyTUulls9DYDbIYSA67UJUj4ZLXs0FtFL9E1nWBkCbo0QArVCqmJGtPZ6HzVtjT6j0OA6KeC2CCFQixJbq5+myHWnHUNX6F+WWj0NgJshhEDtallPWTdSxocq0UsqV59iYQi4HUII1DrnZdL5CfKRT4xp2w2DGgLuhBACdWRgC6UwRe684Bi0XD9dQgwBd0EIgbrTLECsGC5T26gxWfryE7QQcAuEEKhTznf3XTJY/mqzMSXfqODVSQGrEULAAvc2U7YlyyNXHPHZ+tGrLA0BKxFCwBpN/EX2MPlAe7X3Un3hMRaGgGUIIWAZ52XS5cPl9O3m5Dyj3LB6IMCWCCFgsZgmSmGqvFwh+mTrh69wmRSoa4QQsF6Ij/hwoDalixqfrX9whMukQJ0ihIC7SI9U146Ur+w00zcYxbrV0wC2QQgBN9KlobIlSTqEiM3Sd1/kMilQFwgh4F6CfMT7Cdq0KHXgcj1jD5dJgVpHCAF3lB6pbhwt3/nMHLfWuFxh9TSAVyOEgJvq2EDZmixbBIpeS/WdFxy3uk7K9VOghggh4L78NZERp70So45YqaesNt7cf+OV0pwTjt9s5vghUCOEEHB3ae3U/CR5pszxYqH5atG3LVx50vHyDuN/emoWzgZ4AUIIeIC2wcq/R8kx7ZT/KTCf2moIIVaedPy+0Fg2TDbys3o4wMNJqwcAUCV+mpjdR+sfZk7cYGz+wk/TjOXDZQNfq8cCPB8rQsCTjGunzu6j7byknCwWn/O2FYArEELAk6w86XjnoLl7dHkTP9E/R5+Sb1yptHomwMMRQsBjXH9esLm/Y3Oy7Ndc3XXR0XmRPu8Q5+6B6iOEgGdYfsLx0g5j+fCvd8f4qiJziNbAV5kQobxWZA5arh+8zJVSoDoIIeAZwgLFsmHf2x3jq4qPBmmpbdUdqTKptdovW59RyJsaAnesLkK4efPmESNGREdHJyUlHTx4sA4eEfA+PRorP9wj6quKvs0VqYopXdUdqXLvJdEtU19ziqUhcAfqIoS+vr6zZ88uLCycOHHiY489VgePCNhQeD1l4SDtj/eqj200xq81zpZaPRDgIeoihDExMe3btxdCREdHX7x4sQ4eEbCtxNbqvrGyc0PRLbMyY49psjgEbseVB+pLSkpOnDjx3VsCAgJat27t/NjhcDzzzDNTp0514SMC+KFAKWZEa2Paqk9sMt4/bL4Vr8U0UaweCnBfrgzhsWPHXnnlle/ecvfddz///PNCCIfD8atf/SoiIuLRRx914SMCuJXujZS8RPn+ITNxlT4+Qn0pRgv2sXomwC1VKYSlpaX5+fkFBQWXLl36buocDsf777+/atWqZs2aPfnkk126dFmwYMFN72H69OmVlZWzZ892zdQAqkARIj1SHd1anb7d6LRIfyVGTY9kozhwoyr9v6KgoODpp5/etGnTa6+99t3b//KXv7z44ospKSmapsXHxxcXF9/02+fOnTt37tyGDRtOnz79hRdecMHUAKqskZ+YE6/9c4D2apGZuFrnhdmAGyiOW77f5412794dHR1dWfn1CzqZphkREfHmm2+OGDFCCBEXFzdp0qRHHnnkh9945MiRY8eOOT/WNG3AgAE3vf+pU6eeP39+9OjRzk/Dw8Pj4uJ+ZJ6rV68GBwdXcXi4VllZmZRSSl603RrV+8dfaYo39zte2WU+0UmZ1l314+2bqsXhcBQXFwcFBVk9iE2VlJT4+flpWlX/+arq7dd71f9B9uWXX/7nP//p37+/89P+/fvn5+ffNITt27d37hr9caZp7tix49q1a85PY2JioqKifuTrS0tLq/7fAq5FCK1V7X/8j7UTI5qLZ3f4dM80/y9GT2jOa7PdMYfDUVpaWpUfr6gNJSUlhmFU/d9/YGDgbf/HqlEIAwICAgMDnZ82adJk79691b43IYSU8qGHHnrqqaeq+PUOh4Nfyqwiv2H1IDZVk3/8dweJrGEi+7j563y1X3Plj721pv6unc7LORwORVH44WMVVVXvaEVYpfus9ncGBARUVFRcv7JaVlZ2PYoA3Fxia3VvmowIEd0Wc9wQdlf9ELZo0cI0zVOnTjk/PX78eMuWLV00FYBa5zxuuHqE/NdR875l+u6LxBA2Vf0Q1q9ff/DgwfPmzRNCXLp0KTs7Oy0tzXWDAagL3Rspm5Lk4x3VoSv0KfnGVd7dEPZTpRBeuHChffv2I0eONAyjffv2ffr0cd4+c+bM119/fejQodHR0SNHjrx+OwAP4jxuuHesT5khOi3SFx5jBw3spUqbHRo0aJCbm/vt93yzRSI6OvrQoUMFBQXNmjXr1KlTrQwIoE44jxt+8qXjF5uMeYfM2X20NkG8MBtsoUoh1DQtIiLipn8UFBR0/QQFAE93X6iyI1W+sc+MzdJ/0Vl97h7Nl2MC8Hb8GwfwPT6qmNJV3ZIst59zdF+srzvNJhp4OUII4CYigpWcYfK1e9VHPjHSNxjnyqweCKg1hBDALSW2VveOlS3qiS6LOG4Ir0UIAfyYelLMjNVyR8oPj5r9c/Q9l4ghvA0hBHB7UY2UzUly0t3qkOX6lHzjGscN4UUIIYAqcR433DnG51K56LRIX8xxQ3gLQgjgDjQPEPMStPkJ2vMFZuJq/T/XuFIKj0cIAdyx/mHKzjEyvrkavUSfUWhUsDiEJyOEAKrDRxXPRqnbkuXWs47YLH3zGZaG8FSEEED1tQ9RVgyXL8WoP11vpG8wznPcEB6IEAKoqcTW6r6xskU90XlR5dwDnDaEhyGEAFzAedxw9Qj5j8/M/sv0vRw3hOcghABc5p7GyqZE+djd6iCOG8JzEEIArqQq3x43jMrUl59gaQh3RwgBuAoGJggAABGFSURBVF5ogJiXoP39Pu3prUbiav04xw3hxgghgNqSEKbsGiPjm6s9luizdpk6xw3hlgghgFp0/bjh+i/MmCw9/yxLQ7gdQgig1rUPUVYOly/GqBPWcdwQbocQAqgjia3VojGyoR/HDeFeCCGAulPfV2TEaatGyL8fNBOW6fu+ooawHiEEUNd6NFbyk+QDHdT7svUp+UaxbvVAsDdCCMACqiIe76juTvv6uOEKjhvCOoQQgGXCAsW8BO3tftrUrUbiav1EMTmEBQghAIsNCFN2pMqeTZR7MvVZu0yDGqJuEUIA1vPXxIxobWuyXHfajMnSt3DcEHWIEAJwFx1ClFUj5O97quPXGukbjAvlVg8EeyCEANxLYmu1KE029BPdF+vzDnHcELWOEAJwOw18RUaclj1Um73PHJCj7+e4IWoTIQTgpqKbKFuS5P3t1X7Z+rTtRplh9UDwUoQQgPu6ftzwdLHoulhfeZKlIVxPWj0AANyG87jhutOOX2wyIuuLN/pqreopVg8F78GKEIBnGNji6+OGMVl6xh6OG8JlCCEAjxEgxYxoLW+0zDlhxmbpWzluCFcghAA8TGR9ZfUIOf0eNSVXn5xnXK6weiB4OEIIwCONa6fuH+fjr4kui/V5h0yrx4EHI4QAPJXzuOHSIdrre80BOfoBjhuiWgghAM/Ws4mSnyRT2qjxHDdEtRBCAB5PqmJKV7UoTZ4uFt0W66s4bog7QQgBeIkWgcq8BO3PcdoTm4zE1fpJ3t0QVUMIAXiVUa2UvWmyZxOlJ8cNUTWEEIC3uX7ccNkJs9dSfds5YogfQwgBeCfnccMpXdTk1frkPONKpdUDwV0RQgBeSxEiPVLdN9bHXxOdF3HcEDdHCAF4uYZ+IiNOyxqivb7XHMhxQ/wAIQRgCzFNlPwkmdxGvW+ZPqPQKOe4Ib5BCAHYhfO44c4x8ugV0XWxvvoUS0MIQQgB2I3zuOH/9dYm5xnj1xpnSq0eCFYjhADsaHRrZW+a7NxQdM+s5LihzRFCADYVKMWMaG3jaJl93Oy1VN/OcUO7IoQAbO2u+kruSDmli5qcq0/J57ihHRFCAHbnPG64N81HCI4b2hEhBAAhvjlu+OEA7bUic9By/eBlrpTaBSEEgG/Fhyo7UmVSa7VfNscN7YIQAsD3OI8b7kiVey+Jbpn6mlMOIcTeSzdZIJ4sdnxVUefzwdWk1QMAgDsKr6csHKRlHzcf22j0aqoU62JUK+UXnb9dPJwodqTkGn/vp93TWLFwTtQcK0IAuKXE1uq+sbJzQ/HpeXPOAfMve77eR+Os4Bt9qKA3IIQA8GOcxw1zR8h6PuKFHeZvt5knS4Szgvc2o4LegEujAHB73RspmxLlW/vNJ/ONvx/0zR5GBb0HK0IAqBJFiNGtlbsbKFJVxq0xeCUar0EIAaBKnM8Lzu2rHkwqbxUsBi/Xp203Kjh87/kIIQDc3vXdMfc2U3xUkTdaxoeqK447ei7RC8+zNPRshBAAbsMhxEP/Nt7q++3zgj6qWDJEiwhRxrZTRq7Sp203KlkaeixCCAC3oQixcriMbfq93TG+qvhokPZ8tLYj1WfvJUd8tn7gK5aGHokQAsDt+dzsh6WPKhQhwgLFx0Plo3er9y3TZ+3irQ09DyEEgJpShHi8o7otWa46afbL1g/xgt0ehRACgGu0DVbWjpIP3aX2ydZn7TJNaughCCEAuIxzabglSeacMIet1E8UE0MPQAgBwMXahyjrRsrBLdSeS/S5B9hO6u4IIQC4nlTFs1Hq+lFy7gFzxEr9FEtDN0YIAaC2dGmobEmSCWFqdBZLQ/dFCAGgFjmXhrkj5Bv7zPFrjfNlVg+EHyCEAFDrujdStiTLiBDRPbNyyecsDd0LIQSAuuCviZmx2uLBctp2c/xa42K51QPhG4QQAOpOXDNlZ6qMCBHdFuvZx1kaugVCCAB1KkCKmbHavwZq/2+Lmb7BuFpp9UC2RwgBwALxoUphqgyQonumvv4LDldYiRACgDVCfMSceO2tvtqDG4zJeUaxbvVAdkUIAcBKw1oqRWlSCBGVqW/8kqWhBQghAFisga+YE69lxGn3rzembTfKDasHshlCCABuYVQrpSBFHrosembpBedZGtYdQggA7qJZgFg8WHshWh21Sp+23ajgeEWdIIQA4F7GtVN3jvHZ/5WIzdJ3XmBpWOsIIQC4ndAAsXSI9lQ3degKfdYu06CGtYkQAoCbSo9Ut6fI1afM+Gz94GViWFsIIQC4rzZBypqR8uG71H7Z+qxdpkkNawEhBAC3pgjxeEd1S7JcfsLsn6MfvkIMXYwQAoAHiAhW1o+SY9uqfT7W5x5gZehKhBAAPIOqiCld1Q2j5dsHzBEr9ZPF1NA1CCEAeJLODZT8JDkgTI1eos89wElDFyCEAOBhpCqejVLXjpRv7jdHr9K/KLF6IA9HCAHAI3VrpGxLljFNlZ5ZlYuPsTSsPkIIAJ7KRxUzorUlg+V/fWqOX2tcKLd6IM9ECAHAs93bTNmRKiNCRPfF+sf/YWl4xwghAHi8AClmxmofDdKe2mqOX2tcYml4JwghAHiJvs2VXWNkWKCIztLXnuZwRVURQgDwHoFSZMRpc+K1Rz4xJucZ1yqtHsgTEEIA8DZDw5WiMVIIEZWp//sLloa3QQgBwAvV9xVz4rXX+2g/22BMzjNKdKsHcmOEEAC81shWStEYWaKL2Cz90/MsDW+OEAKAN2voJ95P0Gb0VBNX6dO2GxUcr/gBQggA3m9cO3XnGJ8DX4mYLH3HBZaG30MIAcAWmgeIrCHa73qow1fqMwoNgxp+gxACgI2Ma6fuTPX59Jyjb7Z+4CtiKAQhBAC7CQsU2cPkI3ep9y3TZ+0yeZNfQggAtqMI8XhHdWuyXHnS7LdMP3TZ1jEkhABgU+2ClXWj5IORap9sPWOPfVeGhBAA7Mu5NPxktFxwxBy+Qj9RbMcaEkIAsLtODZTNiXJgC7XnEn3uAdudNCSEAAAhVfFslLpulJxzwBy1Sj9dYqOlISEEAHyta0NlS5K8L1TtsUSff9guS0NCCAD4lo8qno1Ss4fKl3ea49ca58usHqj2EUIAwI16NVV2pMqIENE9szLrP16+NCSEAICb8NfEzFht8WD5223m+LXGxXKrB6o1hBAAcEtxzZRdXy8N9WXHvXMHDSEEAPyYAClmxmrv3qf9crMxOc+4Vmn1QK5GCAEAtzc4XNmdJoUQ3TP1DV941dKQEAIAqiTER8yJ197oq03cYEzOM4p1qwdyEUIIALgDw1t+vTSMytTzvvSGpSEhBADcmQa+Yk689odYddxafdp2o9yweqCaIYQAgOoY107dNcbns8siJksvPO/BS0NCCACopmYBInOw9ny0OnKVPm27UemZJ+8JIQCgRsa1U3eO8dl3ScRm6bsuet7SkBACAGoqNEB8PFSb2k0dslyftcs0PKqGhBAA4Brpkeq2ZLnqpNkvW//sssfEkBACAFymbbCydpR86C61b7Y+a5dpekINCSEAwJUUIR7vqG5NljknzGEr9ePX3D2GhBAA4HoRwcr6UXJwCzUmS597wK23kxJCAECt0BTxbJS6fpSce8AcsVI/VeymS0NCCACoRV0aKluSZEKYGu2uS0NCCACoXVIVz0apuSPkm/vNcWuNc2VWD/R9hBAAUBe6N1Lyk2T7EBGVWbnkczdaGhJCAEAd8dfEzFgtc7Ccvt0cv9a4WG71QEIIQggAqGO9myk7UmVEiOi2WM8+bv3SkBACAOpagBQzY7V/DdSmbjHTNxhXK60chhACAKwRH6oUpMoAKbpn6uu/sOxwBSEEAFgmxEfMidfe6qs9uMGYnGcU6xbMQAgBABYb1lLZnSaFEN0X6598WddLQ0IIALBefV8xJ177Sx/t/vXGlHyj3Ki7hyaEAAB3MaqVsiNVnioRPbP0T8/X0dKQEAIA3EhTf7FokPZCtDp6lT5tu1FR+8crCCEAwO2Ma6fuHOOz/ysRm6Uv+dx8eeeNPSzVxS83uyaT0gX3AQCAq4UGiKVDtIXHzJ9vMlrXU65UOmbFas4/KjfET9YZae0UX1es5lgRAgDc17h26rZkWd9XvPeZ+XieIYQoN8S4tUZaO+XBSNckjBUhAMCttQlSckfKuQfMJ/ONfRdlAz/HuAjVVRUUrAgBAO5PEWJyR3VLktx7WZGKw4UVFIQQAOARyg3xuwJzVrQRFqhM2+7KY4aEEADg7q4/L/izdubrccqVCuHCFhJCAIBbu2F3jCLEX/tqLmwhIQQAuLVTJY6ftv/eHlFFiNl9tKb+ikvOEXpqCMvKyr744gurp7CvCxcuXLlyxeop7Ov48eOmaf3bmdpTSUnJmTNnrJ7CXiKClZ+2/7pW58+fv3r1qhBCVcRT3VRbnyMsKCh48MEHrZ7Cvl566aUFCxZYPYV99evX76uvvrJ6Cpv65JNPnnjiCaunsK/nnntuyZIlrr1PTw0hAAAuQQgBALZGCAEAtqY4HHX9XsC3smDBgmeeeSYgIKAqX1xZWVlcXNygQYPango3dfXqVSllFf/HgstduHChUaNGiqJYPYgdVVRUlJaW1q9f3+pBbOrq1as+Pj7+/v5V/Pply5Z16tTpx7/GjUJomubnn39u9RQAAO/RsmVLX1/fH/8aNwohAAB1j+cIAQC2RggBALZGCAEAtkYIAQC25nnvUG8Yxv79+3fu3GkYBq+yVvdOnTqVk5Nz8ODBpk2bPvDAA61atbJ6IhtxOByLFy/etWtXcXFx586d77///sDAQKuHsqPCwsJPP/10woQJISEhVs9iI/n5+bt3777+6aOPPqppmkvu2fNWhDk5OcOHD//rX/86depUq2exo4ceeigvL69ly5bHjh3r3Lnz3r17rZ7IRkzTnD9/fmBgYNu2bd9///2EhATDcOXbk6Iqrly5MnHixMmTJ587d87qWexl4cKF77333tFvuPDIg+cdnzBNU1XVzZs3JyYmXrhwwepxbKesrOz6Uda0tLQOHTrMmjXL2pHsqaSkJCQkZM+ePR07drR6Fnt5/PHH77nnnl/+8peHDx9u37691ePYyNSpU4OCgn7/+9+7/J49b0Woqp43szf57gs6lJWVBQUFWTiMnW3evLl+/frh4eFWD2Iv69evP3To0COPPGL1IDa1Y8eOWbNmffjhh2VlZS68W6KCalq+fPnWrVsnTZpk9SC2k5aWFhoaOmbMmI8++ig4ONjqcWykpKTkV7/61ZtvvsmL21kiLCwsPDz8ypUrGRkZUVFRly5dctU9e95mGbiDrVu3Pvjgg/Pnzw8NDbV6Ftt59913r1y5smjRogkTJhQVFYWFhVk9kV1Mnz594sSJHTt2LC8vt3oWO3rmmWecH5im2a9fv9mzZ//ud79zyT2zIsQd27FjR3Jy8rvvvjt8+HCrZ7Gj4ODg8PDwKVOmtGrVKjc31+pxbGTBggUffvhhTExMnz59hBApKSk5OTlWD2VHqqrGxcUdPXrUVXfIihB3pqioaOTIkXPmzBk1apTVs9hOWVmZn5+f87rcuXPnPv/88zZt2lg9lI2sXbtW13UhREVFRZ8+fV566aXevXtbPZSNlJaWOt/xpqysLDc3Nz093VX37Hm7Rk+ePJmSknLt2rUjR45ERUW1a9du4cKFVg9lIzExMYcPH+7QoYPz0yFDhvzhD3+wdiT7yM3N/fnPf96zZ09FUdavX5+UlPT222/zfFXdKy8v9/f3Z9doHWvXrl2XLl0aNGiwcePGdu3arVixwlXvBOd5ISwvL9+zZ8/1T/39/bt06WLhPHazb9++0tLS6582bNgwIiLCwnnsZv/+/fv27VMUpWvXrnfddZfV49iUw+EoLCzs2rWrn5+f1bPYyMmTJwsLC0tKSjp06BATE+PCe/a8EAIA4EJslgEA2BohBADYGiEEANgaIQQA2BohBADYGiEEANgaIQQA2BohBADYGiEEANgaIQQA2BohBADY2v8HJT3g4vyOUtIAAAAASUVORK5CYII=", "text/html": [ "\n", "\n", "\n", " \n", " \n", " \n", "\n", "\n", "\n", " \n", " \n", " \n", "\n", "\n", "\n", " \n", " \n", " \n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n" ], "image/svg+xml": [ "\n", "\n", "\n", " \n", " \n", " \n", "\n", "\n", "\n", " \n", " \n", " \n", "\n", "\n", "\n", " \n", " \n", " \n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n" ] }, "metadata": {}, "execution_count": 5 } ], "cell_type": "code", "source": [ "p" ], "metadata": {}, "execution_count": 5 }, { "cell_type": "markdown", "source": [ "The `info` object passed to the callback contains not just the densities\n", "but also the complete Bloch wave (in `ψ`), the `occupation`, band `eigenvalues`\n", "and so on.\n", "See [`src/scf/self_consistent_field.jl`](https://dftk.org/blob/master/src/scf/self_consistent_field.jl#L101)\n", "for all currently available keys.\n", "\n", "!!! tip \"Debugging with callbacks\"\n", " Very handy for debugging SCF algorithms is to employ callbacks\n", " with an `@infiltrate` from [Infiltrator.jl](https://github.com/JuliaDebug/Infiltrator.jl)\n", " to interactively monitor what is happening each SCF step." ], "metadata": {} } ], "nbformat_minor": 3, "metadata": { "language_info": { "file_extension": ".jl", "mimetype": "application/julia", "name": "julia", "version": "1.7.2" }, "kernelspec": { "name": "julia-1.7", "display_name": "Julia 1.7.2", "language": "julia" } }, "nbformat": 4 }