{ "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", "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 Eₙ-Eₙ₋₁ ρout-ρin α Diag\n", "--- --------------- --------- -------- ---- ----\n", " 1 -7.844783373684 NaN 1.98e-01 0.80 4.0\n", " 2 -7.850314129179 -5.53e-03 2.93e-02 0.80 1.0\n", " 3 -7.850646691099 -3.33e-04 3.02e-03 0.80 3.0\n", " 4 -7.850647500382 -8.09e-07 4.69e-04 0.80 2.5\n", " 5 -7.850647510664 -1.03e-08 1.74e-05 0.80 1.5\n", " 6 -7.850647511581 -9.17e-10 6.05e-06 0.80 4.2\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+gvaeTAAAgAElEQVR4nO3deUATZ94H8GeOhBBu5AZRRFCOeoAHKirI5QEKCh49cHvvduuy7vZwt7Wrbbet7dqufXtp221Lra31QuVQwKOC4okn1HoroqIoIkeCmZm8f8TirQkJmUzy/fyVjMnMj2mab55nnucZSqvVEgAAAFtFi10AAACAmBCEAABg0xCEAABg0xCEAABg0xCEAABg0xCEAABg0xCEAABg0xCEAABg0xCEAABg0xCEAABg08wUhIWFhQkJCUlJSTt27Ljfa65fv37w4EH998nzvClKsyE4Y4YSBEHsEiQGZ8xQOGOGEgTB5CuDUmZYa/Ty5cuDBw/esWNHU1NTUlJSVVWVXC6/+2WnTp2Kj48/efKknrttampycnIyaaVWDmfMUM3NzY6OjmJXISUtLS1KpZKiKLELkYzW1laFQkHT6JzTl0qlksvlDMOYcJ/mOPvFxcVJSUldunTp3r17WFjYAxqFAAAAZmaOIDx//ryPj4/usa+v77lz58xwUAAAAH2YIwgdHBza2tp0j9VqtYODgxkOCgAAoA+jgvDKlSv5+flvvfXWF198cet2tVo9b968adOmvf32283NzeHh4fv379f904EDB8LDw405KAAAgAmxxrz5hx9+WLp0qe7C+B//+Mf27U8//XRdXd0LL7zw3XffTZkyJT8/X6VSzZ49u6GhITw8vEePHh0+YpOGyGiiuOsq6SU18VR0eK8AAGC7TDBq9Msvv8zNzS0rK9M9ramp6dmz59mzZz09PVtaWry9vbdv3x4cHFxSUiKXyxMTE1n23ul76tSpoUOH/ve//23fkpqaqlDclm9rzmgXHRaWJzAK5uYYyK9+0248r10Sh2FXD4FRo4bCqFFDYdSooTBq1FCGjhrV59wa1SK8p127dvXq1cvT05MQ4uDgMGDAgO3bt0dGRo4fP/6h721ra/vpp5/anw4dOtTV1fXWFyR6kAv+THqx8MMwDa9WMQyTe4JZc5ZePEzT2mryP8XaqFQq0445tnoqlQrfUAZRqVSEEASh/lpbWwVBwMdMfyqViuM4/b/KlErlQ0+v6YPwwoULXbp0aX/q4eFx/vx5Pd/r7Oy8cuXKB7/mj48QOzshezvz3WDt0lplwXkhL4VVMHYdr9hmaLVatG8MhTNmEIqi0CI0CE3TaBEahGEYk88jNH0Q2tvbazSa9qdqtVqpVJr2EE+G0lpCEjfadVEI60azdmjkAABAR5n+Z4i/v//p06fbn545cyYgIMDkR+EFouLItjrteweE61iiCAAAOsr0LcKRI0eq1epNmzbFx8fv27fv5MmTo0ePNu0hvv5NyDstbEtp+75GOf8g//0R4dNhTEoAemMAwHxKS0sXLVpk5E54nqdpGp3J+khJSXn66ac7Y89GBeHGjRufffbZpqampqam4ODg0aNHf/rpp3Z2dh999FFWVtbgwYN37do1b948FxcXU5VLCPn6N2HlKWFFIqtpJX97hHGzoz6tFv5Yzke6k8+GMV0d8HkCAHPYuXMnx3HTpk0TuxCbsG3btpKSEksMwpiYmJKSkvan7UvGPP7444mJiVVVVb179/b39zeqwNvlnRbyTgsrE1k7huiuQ+quFxbVaCPcSL+V3Ct9mJf60AzSEAA6X1hYWFZWlthV2ARBEFatWtVJOzcqCJVK5f1mx/v4+LSvL2pCiX70mAD6jtExT4XSE7qRLnbk8Z70i9v4H48Ln8cyQ7wQhgAA8HASG7PrKCP3HCPaxY4QQno6U+tGs28NoKdu5LM38/VqM1cHAADSI7Eg1EdaIH1gIutmR8KXaxYcEoROv98iAABImBUGISHERU4WDGGKx7BLTwgjC7iDVxCGAABwb9YZhDr9ulBbx7PP9qKTi7icCr5J8/C3AAB0qh+PC7lHHzL3uaRW+9EhzI82H2sOQkIIRUh2CF2VKSOEhC/nHvr5AwDoVBO708tPar+7/3dRSa32jT38Ez2t/MvZotjEuXa3IwuGMD/GM/85KCQUcoevoqcUAMRhx5BlCcyK+2ShLgXXJrMeHbqvXLdu3S5fvmxsiYQQQgRBKC0t5XneyP20tbV5eHi035v9nvr27Xv8+HEjD2QMmwhCnVgfqjKdHR9Ijyzg5lTyamP/+wIAdMT9stDIFCSEtLW1GX9nPZ3r168nJSU9OMD0wbLsG2+8cb+777UfSxDE7K4z/RJrloylSU4kndWDmrVTiFzBfTKUGY2F2QDA7HRZmLWBJ0SYHkITU6TgHbRa7dKlS/ft2xccHJydnW1nd+MWPaWlpeXl5V5eXtOmTXNzcyOEbNq0ydPTMzIykhBy9uzZXbt2ZWRk5OXlEUL+97//yeXy8ePH6zkvvLW1deHChfX19UOGDElNTdVtbL+tbF5eXlRU1ObNm48ePZqQkBAXF2eaP9VoNtQibOenpHLjmEWxzMztfFoxd6YZPaUAYG63tgtNnoKEkBkzZixatKhnz56rV69uvx3s22+//eKLL3p5ee3fv79fv35XrlwhhHzzzTcbN27UvaC6uvrdd9/t2BFVKtXw4cMPHz4cGBj4yiuvvPfee4QQjUbz/PPP625JNG/evKlTpx45csTFxWXSpEmbN282/s80CdtqEd5qlB9Vmc7OO8D3X8XNiKD/2Y+R2+KvAgAwpXHrOYPuhyNoSU6FwAlkgCc1bRNn0LGWJbCu8nv/0+nTp5csWXL69GknJ6fp06cHBQWVl5f37dv33//+d2VlZVhYGCEkLS3ts88+e/311++5h/T0dELIU089pf999H7++Wc3N7eFCxcSQmJjY2NiYmbOnHnHa8aMGTN79mxCyJUrV/Ly8iykUWi7QUgIsWfJnCjm8Z70jG38wDzu82HMUG/0lAJAx/39EcagRTz21GvPtQgKlsR4UYl+hv0YV97/+/vgwYNhYWFOTk6EEJlMFh0dfejQIUdHRwcHB10KEkKGDx9eWVlp0BEJIRcvXvz73/9OCGFZ9ptvvrn1nw4cOBATE6N7HBERQVHUqVOnunXrdutr+vXrp3vg7++/bds2Q4/eSWw6CHV6OlNFo9m1Z4RHN/EjfKj5MYyn6XonAMCmjPIz4Md0Sa0277RQlsY6sCSjlAtz1equFxqvpaVFrb65yKRKpVIqlUqlUq1Wa7Va3V2fWltbda09hmHaR4eqVKoH79nR0XHSpEmEEJq+s9SWlpb2x4IgXL9+/e7WpGnvLG8q6A28IS2Qrspk/RxIBBZmA4DOd+t1QXuWrE5i7zenomP2799fXV1NCKmpqdm+ffuwYcOCgoLc3d1XrlxJCFGpVMuWLdP1THbv3l3XNNRqtStWrNC93c7Ozt7eXncR8VZKpTI9PT09Pb39uuOt1qxZ09raSghZsWJF165dTXv3oc6DFuFNDix5byAzLZh+YSu/+Jjw2TBmoCd6SgHA9O4eHXP3OFIj9enT58knn/T29t65c+drr70WHBxMCPnmm2+ys7O/+uqrI0eODBky5LHHHiOEPP3000OGDBk5cmRzc3N4eLju7RRFvfjii9HR0V27dv3qq6/auzQfLCQkZOjQoQEBATt37vzpp5/ubjVaJgThnfq6U+Vp7PdHhQklXFYQ/dYAxlkmdk0AYEXuN0bUtFno4eGRn5+/b9++gIAAPz8/3caEhIRjx479+uuvnp6eXbt21W0MCAg4fPhwdXV1QECAh4dHe5/q+++///bbb7e0tOiuNeojNjb2T3/605EjR8LCwpydnQkhCoWirq5ON4OiuLi4fSrFM88884c//MHIv9FUEIT3oFuYLS2QnlPJhy/n3hlAZ5uo4x4AbNzGc9o39vD5KWwXu3v8qx1Dlo5iMko5liKPGb3KmlwuHzRo0B0b7e3to6Ki7tjo4OAwcOBA3WOZ7OZvf7lcLpffZ2Tqfbi5uQ0ePPjWLV5eXroHtwaqnZ1d+9RG0eH7/b7c7MiCIUxeEvN/VcKoAu5XLMwGAEYLdib3S0Ede5asSmSjPDpyXaZr16660ShhYWFTpkzpcJEdk5KSMmTIkA680d/f/9b0NT+0CB9igAdVMZ79tFoYvpZ7pjc9J4pRWOKgJwCQhm6OD084e5aEuXYkCHft2qV70KdPnz59+nRgD8bIyMjo2BtLS0tNW4mh0CJ8ON3CbAcnyc61kMgVXGENmoYAANYDQagvXyXJjWO+HM68tINPK+ZOY2E2AACrgCA0TLwvtX8iG+tNR63i5lTyBq2lBAAAFghBaDAZTV7tS++cwO68pO2zgttwDk1DAAAJw2CZDgp2pgpT2LVnhGfK+OHe1H8GM172YtcEAGYUFBT00ksvrV+/3pidtC94Bg/W0NDQvpCpySEIjZIWSCf40e8f4B9ZqflnX+bFCJrBRxrANkybNq1Xr15G3ghXrVbL5XKprMAiroCAgE7aM4LQWEqWzIlipgXTf97K5x4TPh/GDMLCbAC24e6Z6YZqbW1VKBQIQnHh7JtGLxeqZCybE0GPL+aeL+cbr4tdEAAA6AdBaDK6hdl+zZQpGBKxgss9ijtYAABIAILQxHQLs61OYj6pFuILuGoszAYAYNkQhJ0i2oPaPp59NJgesZbLqeCbNWIXBAAA94Eg7Cw0RZ7rTR+cJGtoI31XcgVYmA0AwCIhCDuXbmG2r0cwr2BhNgAAi4QgNIc4X2rfRDbRjx6Yh4XZAAAsC4LQTGQ0yYmkt09gd1/SPrKCK61F0xAAwCIgCM2qhxOVn8L+ZzD9TBk/eQNfpxK7IAAAm4cgFEFaIF2dyYa7kT4rNQsOCTwahwAA4kEQikO3MFtZKptfIwxaze24iDAEABAHglBMoS5UyRh2Vl86vYTL3sxfbhO7IAAA24MgFF9WEP1rlszNjvTBwmwAAGaHILQIrnKyYAizNpn57FchLp+rakAaAgCYCYLQgkR5UNvS2Kd70QmFWJgNAMBMEISWhaZIdgi9b6KsoY2ELeeWncTcewCAzoUgtEQ+9iQ3jlkcx8ytFNKKuVNN6CkFAOgsCELLNdKX2pvBJvrRg1Zzcyr5Nl7sggAArBGC0KLpFmbbm8FWNZBHVnLFWJgNAMDUEIQS4O9ALUtg5g+mny/n04q5sy2IQwAAk0EQSkZaIF09iY32oKLzOCzMBgBgKghCKbFnyZwopjyVLagRBuRx27EwGwCA0RCE0hPiQhWPYd+Mpidv4LM38/VqsQsCAJAyBKFUpQXSByaxbnYkfLlm0WEszAYA0EEIQgnTLcy2fgz7zRFhZD53CAuzAQAYDkEoef27UFvT2Gd60UmFXE4F34SF2QAADIEgtAa6hdmqMmVqnoQt53KPYmE2AAB9IQith7sdWRjL/BjPfHBASF3PncTCbAAAekAQWpvhPlRlBpvkTw9ezc2p5JefFN7ce2cD8bdG7WObsGIbAAAhCEKr1L4w24lr5B+7hB0XtXMrb2bhkUbt1I38P/vhPz0AACGEsGIXAJ3F34HKjWPyz2j/UsFXNQiN17VzI8iRRu2UjfziOCbCjRK7QAAAi4AgtHKpgdQoP/adffz8g8LOOnmzwP+AFAQAuAX6x6yfkiVvD2CWJzJ7rtD1auIqF7sgAABLgiC0CUcata/vFjYltXVzIL2W4cb3AAA3malrtLy8vKyszN7e/q9//at5jgjt2q8LBrLa8vFsZin/x3Jh7WntF7GMEl3jAGDzzNQirKiooCjqf//7n3kOB+3uGB1DEbI8kZnQjdp/RTswjzt4BXMNAcDWmSkIX3755UcffdQ8x4Jb8VqyJP620TEUIV8NZ+ZG06/2pUcVcgsOoZsUAGwausasXJjrPQaI0hRJ70YTQgZ7UdM28lsuaL8czrjbmb04AAALYOIgTE5OvnLlyq1bCgoKvL29TXsUMJVeLlTFePZflXzUKm5xHBPrg2kVAGBz9A1CjuOqqqqOHj0aExMTEBDQvr2hoWHFihUqlSotLa179+7FxcWdUyd0FjuGvDeQSfDTTtvEP92Lmt2fYZCGAGBL9L1GGBwcPGHChOnTp1dUVLRvvHLlSlRUVGlp6bFjx/r3779///77vb2ioqKgoKCxsXHZsmWnTp0ysmgwuSR/anc6u+OiNrGQq23BCBoAsCH6tgj379/v6ur6yCOP3Lrx66+/Dg0N/emnnwghjo6O8+bNW7JkyT3ffvHixWvXrr3wwgsnTpzo3bv3/Y7C8/yJEyfan3br1o1hGD0rBCN525PC0ezHh4SBq7mFsUxaIOaYAoBN0DcIXV1d795YXFycnp6ue5yWljZ27Nj7vX3ChAkTJkx46FEaGxsTEhLan27YsMHT0/N+L25ubn7oDuFW+pyxp7qRfk70k9tkS34TFgzk7Bmbbh22tLRotTZ9BgzV2trK8zxFoXtdXyqVSqPR0DR+d+pLpVLJ5XL920hKpfKhLzZqsExtba2vr6/usZ+fX0NDg0qlsre37/AO3d3dT548qf/rnZycOnws26TPGRvpRPb5kT+V86M2sD/GM4+42+6XGkVRjo6OYlchJTRNK5VKBKH+GIZRKBQIQv2xLGtQEOrDqLNPUVT772X8cLYmzjLyQzzzSh9MNAQA62dUEPr6+tbV1ekeX7hwwc3NzZjmIFia7BC6PI397qgwqZRvaBO7GgCAzmFUECYnJ69du1b3OD8/Pzk52RQlgQXRTTQMcCD9V3HlF9DoBwArpG8QvvXWW5MnT66pqfnoo48mT57866+/EkKefvrpX3/99bHHHps5c+Ynn3wya9asziwVxGHHkAVDmEXDmamb+DmVPI80BADrou9gmZEjR/bu3TsrK0v31MPDgxDSpUuXysrK5cuXq9XqysrKoKCgzioTxJbsT+1JZ6f/wiUVct/HMf4OGA0BAFZC3yAcMWLEPbe7u7s/99xzpqsHLJe3PSnCREMAsDr4LgMDUITkRNJ5SezftgvPl/MqTuyCAACMhiAEgw3ypPZksM0aMnA1d6gB1wwBQNoQhNAR7RMN4wsw0RAApA1BCB2XHUKXpbLfHhUyN2CiIQBIFYIQjNLbldo+nvVXkv6ruK116CYFAOlBEIKxdBMNPx5CTyrlMNEQACQHQQimMb4bvX+irKJOm4Q7GgKApCAIwWS87cm6MeyEbvTA1Vz+GWQhAEgDghBMqX2i4V+38zkVfBsvdkEAAA+DIATTG+RJVWawl9QkOg8TDQHA0iEIoVM4y8gSTDQEAClAEEInwkRDALB8CELoXL1dqYrxrL+SROVhoiEAWCIEIXQ6BUMWDGH+G3NjoqGANAQAS4IgBDOZ0D7RsIg714owBABLgSAE89FNNBwfSA/M4zHREAAsBIIQzEo30XBVEoOJhgBgIRCEIIJBntSedPaimgxdyx1pRNMQAMSEIARxuMjJj/FMTgQ9bC0mGgKAmBCEICbdRMNvjghZG/ir18WuBgBsEoIQRNbbldo+gfVTkv6ruG2YaAgAZocgBPG1TzSciImGAGB2CEKwFBO60bvT2U3nMNEQAMwKQQgWJMCB2jTuxkTDghpkIQCYA4IQLAtNkZxIemUik1OBiYYAYA4IQrBEg70w0RAAzARBCBaqfaLhiHxu0WFMNASAzoIgBIuWHUJvGsd+Vi1MxkRDAOgcCEKwdGGu1PYJrC8mGgJA50AQggToJhp+hImGANAJEIQgGend6F3p7KZz2mRMNAQA00EQgpR0daA2jmOT/OmBeXwhJhoCgCkgCEFiGIq82pdekcjM2MbnVPDXMZ4UAIyDIARJivGiKjPYOhUZuoY7iomGAGAEBCFIlYuc/DSK+UsEPRwTDQHACAhCkDbdRMNPMdEQADoKQQiSF+ZK7ZjA+ipJ1Cqu4iK6SQHAMAhCsAa6iYYfxtAZJZhoCACGQRCC9bh1ouH5VrGrAQCJQBCCVdFNNIz1oaJWaTDREAD0gSAEa8NQZE4UsyqJxURDANAHghCsU/tEw2GYaAgAD4QgBKulm2g4I4Iens99fwwNQwC4NwQhWLnsEHrjOPY/BzDREADuDUEI1i/8lomG2zHREABuhyAEm6CbaDh/MD0BEw0B4HYIQrAhGd3p3ensxnPalHWYaAgANyAIwbZ0daA2jWOHeVPReZoiTDQEAAQh2CDdRMMf49nnyzHREAAQhGCrRvpSByexF1Rk2Bru2DU0DQFsF4IQbJeLnCwdxcyIoGPXcosx0RDAViEIwdbpJhq+f0DI3sw3a8SuBgDMDkEIQMJdqZ0TWDc70mcl978jwtozd7YO1Tx5dz/mXABYJ1bsAgAsgm6iYZyv8Hw5b88SXiDp3W/8TGzjyeQN/LhAihK3RADoHGgRAtykm2jop6SeLOO/PSIQQtp4krWBHxdIPd8b/7MAWCe0CAFuE+hIlaWys/fwz5fzVb2Z35qRggBWDkEIcCeWJu8OZAZ60lM3kt4u2swgRuyKAKAT4XcuwD208eTbI8Jrj/CX27Q9ftIsOCTwGCoDYKUQhAB3ar8u+PcwrjpTFuJCfXOUH5jHbatDGAJYITMF4b59+z755JOFCxdevHjRPEcE6Jg7Rse4yMmGsayzjBrblXp0E5+9ma9TiV0iAJiUOYJQrVa//vrrdnZ2LS0tgwYNunTpkhkOCtAxOy9pJ3a/bXSMi5ysTWZVPKnKZHs4kz4rNQsOCRwWogGwFuYYLKNQKPLz83WPt23btnv37jFjxpjhuAAdMNyHGu5z54xBFzmZP5ghhMyJYh4Npv9SwX9zRPhkKBN71ysBQHLMeo3w8uXL+/fv79+/vzkPCmBaoS7UutHsWwPoJ37hJ2/ga1pw4RBA2kzcIpw+ffqpU6du3fL111/37NmTEKJSqaZMmTJ37lwfHx/THhTA/NIC6QQ/+v0DfL+V3IwI+p/9GDlGngFIkwFByPN8TU2Nu7u7s7PzrdsPHDjQ2Ng4cOBAhULx+eefC8JtF0+USiUhpK2tbdKkSVOmTHn00UdNUjeA6JQsmRPFPN6Tzqng+6zgPh7KJPujpxRAevT9ETtt2jRnZ+fg4OAlS5a0b+R5fuLEiRMnTnzjjTdCQ0OPHj2qVCodb0fTNMdxWVlZ4eHhCQkJJ06caGpq6py/BUAEPZ2pghT2g8H0H8v5tGLuTDN6SgEkRt8gnDFjxvHjx0eMGHHrxoKCggMHDuzfv3/Tpk2ZmZn/+te/7vne1tZWhUJx5syZWbNmzZo1a8+ePcZWDWBh0gLpqklstAfVfxU3p5Jv48UuCAD0pm/X6NChQ+/euGLFikmTJjk4OBBCpk+fHhMTIwgCTd8Zrs7Ozj///LM+R7l8+bKbm1v70927d3t5ed3vxc3NzXqVDr/DGTNUS0uLVmtAC+/vIWSiL/XqXjbiKP9+fy7J1+bysLW1led5Cjfq0JtKpdJoNHd/bcL9qFQquVzOMPoufKhUKh/6YqMGy9TU1AwePFj3uHv37mq1+uLFi8aMhXF3d9+7d2/7UxcXlwd/PpycnDp8LNuEM2YQiqIcHR0NeksfJ1LkS0prtTMq6G9Okf8bwnR3sqFUoGlaqVQiCPXHMIxCoUAQ6o9lWYOCUB9GnX21Wi2Xy3WPdQ9UKqNW3aAoyu0W+HCARCX6U/snsol+9KDV3JxKXm1zLUMAKTEqaXx8fK5cuaJ7XF9fr9tigqIApE9Ok5xIujKDPXGNPLKCK6jBIBoAC2VUEA4ePHjLli26x2VlZX379rW3tzdFVQBWIsCByo1jFsYyr+zg04q5E02IQwCLo+81wtWrVx8+fPjMmTPFxcWNjY0ZGRmhoaFPPvnkBx98MHfu3N69e7/yyivvvfdep9YKIFGj/Kh9E9nPqoUha7gnQ+nX+zGOMrFrAoDf6dsibGlpaWhoyMrKCg0NbWho0Gg0hBAvL6+tW7deunQpPz9/wYIFjz/+eGeWCiBhMprkRNJ7M9hzLSR8OZd7FIt2A1gKyqDR4Z3q1KlT8fHxJ0+e1PP1TU1NGANpEJwxQzU3Nxs6alQfm89rZ2zjve3J/w1lwlytaoBlS0sLRo0aRDfNGgMD9Wfo9Al94OwDmFucL1WZwaYF0sPXcjkVfJNG7IIAbBuCEEAEup7Sg5NkDW0kDD2lAKJCEAKIxldJcuOYJfHM/IPCqAKuqsFSrlMA2BQEIYDIRvhQe9LZqcF0YiGXU8FfQ08pgHkhCAHEx9Lkud70oUwZISRiOZd7VEDbEMBsEIQAlqKLHVkwhFmVxHxaLcTlcwevIA0BzAFBCGBZBnhQFePZp3vRKeu4nAq+8brYBQFYOwQhgMWhKZIdQh+aJCOERKzgFh0W0FUK0HkQhAAWyt2OLBjCrElivj0ixKzhdl5CGAJ0CgQhgEWL8qC2jmdfDKcnFHPZm/l6tdgFAVgdBCGApaMIyQ6hqzNlbnYkfLlmwSH0lAKYEoIQQBrc7MiCIUzxGHbZSWHQam77RYQhgGkgCAGkpF8XqiyN/UsEPbGUy97MX0JPKYDREIQAEnNrT2nEcs2CQwKPxiGAERCEAJLkKicLhjBbUtmCGmFgHretDmEI0EEIQgAJ6+1KFY9h50bTj27iszfzdSqxCwKQIAQhgOSlBdJVmWwPZ9JnpWbBIYHDPZ0ADIEgBLAGDiyZE8WUpbJFZ4UBeVz5BfSUAugLQQhgPUJdqHWj2bcG0I9v5tOKuZoWxCHAwyEIAaxNWiBdnclGe1D9VnJzKvnr6CkFeCAEIYAVUrJkThSzYwK765K2zwpu/Vk0DQHuC0EIYLV6OlMFKewHg+k/beXTirkzzYhDgHtAEAJYubRAumoSG+1B9V/Fzank23ixCwKwMAhCAOtnz5I5UczOCeyeeu0jK7miGjQNAW5CEALYimBnam0y+9lQ5m87+LRi7lQT4hCAEAQhgK1J9Kf2T2QT/ehBq7k5lbwaPaVg8xCEADZHTpOcSLoygz1xjUSu4ArQUwq2DUEIYKMCHKjcOGZRLPPKDj6tmDuBnlKwVQhCAJs2yo/ap+spzeNyKvhmjdgFAZgdghDA1slokhNJH5jENsnBY64AABphSURBVLSR8OVc7lEsRQO2BUEIAIQQ4qekcuOY3DjmgwNCYiH361X0lIKtQBACwE1xvtTeDDYtkB6+lsup4JvQUwo2AEEIALdhaZITSR+cJGtoI2HoKQUbgCAEgHvwVZLcOGZJPDP/oDCqgKtqQE8pWC0EIQDc1wgfak86OzWYTizkcir4a+gpBWuEIASAB2Fp8lxv+lCmjJAbY0rRNgQrgyAEgIfrYkcWDGHykpjPfhXi8rmDV5CGYD0QhACgrwEe1LY09uledHIRl1PBN14XuyAAU0AQAoABaIpkh9BVmTJCSMQKbtFhYWud9t/77hxZermNPLmFRy8qSAKCEAAM5m5HFgxh1iQx3x4RZm7nK+q0b+69mYUNbWRCMZcZRNOUiDUC6AtBCAAdFOVBbR3PvhhO76kXco8Is3YJhJCGNpJWzP2jHzOuK2IQpIEVuwAAkDCKkOwQOjWQfm0X//Ehfu8ltoXnkYIgLWgRAoCx3O3I57FMfgq7pY6+qCJjkYIgKQhCADCBhjbyxh7+qxjumoYMWcOJXQ6AARCEAGCs9uuC6YFC9STmWCOJL0QWgmQgCAHAKHeMjnGzI1WZ7L56kryOF7s0AL0gCAHAKHUq7Zyo20bHeNuTAxOZvfXaz6tx5wqQAIwaBQCj9HalerveubGrI7UnnRmRz9ux5KlQ/OAGi4YgBIBOEehIFY9h4gt4JxnJCkIWguVCEAJAZwl1oYpGM8lFnCNLjcGcCrBU+JkGAJ2ojzuVl8T+YQu35QIWHgULhSAEgM4V40X9GM9mbeB21yMLwRIhCAGg043yo76MZdLWc1UNyEKwOLhGCADmML4brRFIchG/cRzTywXXC8GCIAgBwEwmBdHNHEku4rekMt0ckYVgKRCEAGA+00Poa9dJUhG/JZX1sRe7GgBCCK4RAoCZzYigH+9JJxdxV9rELgWAEIIgBADze6M/PbYrNWYd16QRuxQABCEAiOLdgUy0BzVmHdeC21SA2MwUhOfOndu8efPu3bsFAYvwAgChCPl0GNPLhcoo4dpwmwoQlTmCUK1WP/XUU2vXrv3ggw/i4uI0GvSGAAChCFk0nHGVU9M28Rx+IYN4zDFqVKFQrFu3Tvd44MCBx44dCwsLM8NxAcDCMRRZHM9klHBPl/HfjGBoTKkAMZipa1QQhNLS0k8//dTFxSU4ONg8BwUAyyenyfIE9nSz9i8V6CEFcZiyRchx3CuvvHLHxg8//JAQwvN8aWnp2bNnPT09TXhEALAC9ixZk8wmFnJ/285/GMOIXQ7YHBMEoVqtJoQoFAqGYdLT0+/5GplM9t577xFCpkyZUlRUNGHCBOOPCwBWw1lGisew8QXcu/uFf/TFaHYwK72CkOf55557bs+ePWfOnCkpKYmOjtZt12g0zzzzTF5eHiFkypQpn3/++YgRI+5++7Vr1+RyuUKhUKvVJ06c8Pb2NuEfAADWwVVO1o9mR+RzdjT52yPIQjAffT9tYWFh//3vfzUaDcfdnPXz1VdfHTx48Ny5c2fPnq2oqFi8ePE931tTUxMbGxsdHR0TE5OVlRUTE2OCwgHA6njZk5KxzCfVwle/YRQpmI9eLUKGYV566SVCCE3fFpy5ubkvvPCCg4MDIeT555/Pzc2dPn363W+PiIjYvXu3PgdqbGwcNWpU+9PvvvvOzc3tfi9uaWmhKAwyMwDOmKFaWlrELkFiWltbBUEw5mPmRkjeSGrMRjnLqTO7WX8ctra2chx3x1crPIBKpZLL5Qyj77VkpVL50NNr1DXC48ePt0+E6N279/Hjx43ZGyFEqVT+4x//aH/q4eFhZ2d3vxfzPK9UKo08ok3BGTOUIAg4YwbRarVKpdLI31uRSlKYrB29XtbFkRnX1VSlWS6FQoEg1B9FUQYFoT7n1qggbGxsbP+acHR0vHr1qjF7I4TIZLKkpCQ9X0zTND49BsEZMxTOmKF0Z8z4joe+HmR1Mp1WzC0dxcb5WnM3Bv07sQuRjM44Y0bty9PTs7GxUff46tWrmBoBAKYyyJNamchO2ciVX8BN7aFzGRWEERERlZWVusd79uyJjIw0RUkAAIQQMsybWhzHTtrAVdYjC6ET6ds1um3bNt1F3V27djU1NQ0bNsze3v7555/PycmJi4vjOO6TTz759ttvO7NUALA5Sf7UF8OY1GJuw1g2zNWa+0hBRPoG4aJFi2pra4cOHbp69erVq1cvXrzY3t5+4sSJtbW1Tz31FEVRc+fOTUlJ6dRaAcAGZXSnmzmSUsT/ksoEOSELwfQordZS+hxOnToVHx9/8uRJPV/f1NTk5OTUqSVZGZwxQzU3Nzs6OopdhZS0tLQYP2r0nj6rFj48JGxJZfyUVpWFra2tGDVqEEOnT+gDZx8AJOCFcPrP4fSoAr5OJXYpYHUQhAAgDTMj6YndqdHruIY2sUsB64IgBADJeGcgk+BHjVvPNeP23mA6CEIAkJIPBjOR7lRGKafG7QvBRBCEACAlFCGfD2M8FNTUjTxn/WuRgjkgCAFAYhiK5I5keK32yS28YCnD3kHCEIQAID0ymvw8iq1t0T5TxiMKwUgIQgCQJHuW5KewR69pZ27H1UIwCoIQAKRKyZL8ZLbsgvbNvbhaCB2HIAQACXORk3Wj2Z+OCx8cQBZCByEIAUDaPBWkeAzzxa/CwsPIQugIo27MCwBgCQIcqOIxTFwB7yQjjwbj9z0YBkEIANYg2JkqHcuMKuDlNMkMQhaCARCEAGAlerlQhaOZlCLOgaXGdLWqm1RAp8LvJgCwHn3dqVVJ7B+2cGUXML0Q9IUgBACrMsSL+iGOzdzA7alHFoJeEIQAYG0S/alFsUxaMVd9FVkID4drhABghSZ0ozUCGbOO3zSO6eGE64XwIAhCALBOmUF0k4aMKuB/SWW6OSIL4b4QhABgtZ4MpZs0JKmI35LK+tiLXQ1YKlwjBABr9pcIeloPKqWIu9ImdilgqRCEAGDl5kYzKQHU2PVck0bsUsAiIQgBwPrNG8T0c6fSSzg1btkEd0EQAoD1owj5bBjjY09NKObakIVwOwQhANgEmiK5cYyznHp0E8/hNhVwCwQhANgKhiI/xDMqXvtMGS9gqj38DkEIADZETpPlCezJJm1OBXpI4QYEIQDYFiVL1qaw2y9qZ+9BFgIhCEIAsEHOMlI0ml11SvveflwtBAQhANgkDwXZMJb99ojw0SFkoa1DEAKAjfK2J8VjmI+rhK9/QxbaNKw1CgC2K9CRKh7NxBfyTjIyuQcaBjYKQQgANi3EhSpMYVKKOAcZNa4rblJhi/ALCABsXR93Ki+JfWoLt/k8ZhfaIgQhAAAZ7EX9NIqdspHbdQlZaHMQhAAAhBAS70t9H8emFXN7LyMLbQuCEADghmR/6rNhzLj13OGryEIbgsEyAAA3TexON2lIyjr+l3FMdyeMnbEJCEIAgNtMD6GbNCSpiN+Syvoqxa4GOh+6RgEA7vRiOP2HUDq5iLvcJnYp0PkQhAAA9/BaPzo1kEos5K5eF7sU6GQIQgCAe3t3IBPvS41bz7VwYpcCnQlBCABwX/NjmDBXKqOEa8Mtm6wXghAA4L4oQhbGMu521NSNPIelua0UghAA4EEYinwfx2gE7ZNbeAHTC60RghAA4CFkNFmWwJ5t0c6oQA+pFUIQAgA8nD1LViezuy5pZ25HFlobBCEAgF6cZaR4DPvLee2/9+FqoVVBEAIA6MtVTtaNZhcfE/5zEFloPRCEAAAG8LInJWOYz6qFLw8jC60E1hoFADBMgANVPIaJL+CdZGRqMJoTkocgBAAwWE9nqmg0k1TIOchIWiCyUNrw3w8AoCMi3aii0eyzZfz6s5hdKG0IQgCADurXhVqZyGb/wpVdQBZKGIIQAKDjhnpTP8SxmRu4ynpkoVQhCAEAjJLoTy2MZVKLueqryEJJwmAZAABjpXejmzVkdBH/SyoT5ESJXQ4YBkEIAGACj/ekrwskvoDfksoEOiILpQRdowAApvFUKP3XSDqpiK9TiV0KGAJBCABgMn+NpCf3oFKKuIY2sUsBvZkvCA8fPuzj4/PTTz+Z7YgAAOb3VjST5E+NXc81a8QuBfRjpiAUBOGll16KjY1ta8PPJACwcu8PZvq4U+klnBq3bJICMwXhggULJk2a5OfnZ57DAQCIiCLk82GMlz01ZSOvwdLcFs+Uo0YvXLhQWlp66xYvL6/k5OSTJ09u3LhxzZo1e/fuNeHhAAAsFk2R70Yyk0r5aZv4paMYBsNILZgpW4QajebS7a5evUoImTVrVlhY2Jdffnno0KHNmzcfOXLEhAcFALBMMposT2RaNNpnynjMtLdkerUIL168+O233+7Zs6exsXHdunXt21Uq1auvvlpcXOzp6fnmm2/Gx8fPnDnz7rc/8cQT586dM1nJAAASIafJikR29DrurxX8giGM2OXAvekVhOfPnz927FhISMi8efNu3T579uxDhw4VFhbu3LkzPT396NGjXl5ed789NTVV9+DQoUPR0dGhoaHG1w0AIAlKlqxNZhMKuTmV/JwoZKElorRafZvsBw8ejIqK0mhujAjWaDTe3t6FhYUxMTGEkLFjx44aNeqll156wB5aW1tZlpXL5ff811OnTg0aNGj69OntW1599VUnJ6f77a2pqekB/wp3wxkzVHNzs6Ojo9hVSElLS4tSqaQoXBC7U30blbSeejxY+/eI275yW1tbFQoFTWNKt75UKpVcLmcYfX9SyOXyh34gOz5Ypra2tqGhITo6Wvc0Ojr60KFDD36LUql88AtomnZzc2t/yjDMAz4fNE3j02MQnDFD4YwZSnfGEIR387InRckkvohykVPP9bq5nf6deKVJTGecsY4HYX19vVKplMlkuqeurq6VlZVGVmNvb//Pf/5TzxfLZLL2o4M+cMYMhTNmKN0ZQxDeUzcXsjlVOyKft5fTT4Xe+B7XnTEEof44jpPJZPq3CPXR8bPv6uqqVqt5/sZ80aamplsbcwAAcIdAR6p4DPP6buHrI/eYXXi+1fwVASHGBKGfn59MJmufC/Hbb7/16NHDRFUBAFinUBfq/4bSfyrjfzh228XCdWe1kzdyHGbfi0GvINRqtQ0NDdeuXSOENDQ0NDY2EkKUSmVmZub8+fO1Wu3hw4cLCwsff/zxzi0WAED6JnWnPx3G/GELt/LUjdxbf1b7ZiWfl8Sy6CIVg17XCC9fvqyb8+Dk5BQcHOzr61tVVUUI+eCDDyZPnuzh4aHVat99913MiwAA0MezvWk1T6Zs5JeNpBRy8uZefm0K28VO7LJslV5B6OHhceXKlbu3+/r6lpWVtbS0KBQK0166BACwbjMi6FZOm7WZDXTkt6TKkIIiMkE73MHBASkIAGCovl3oro7aS2oSvlwzZSO/4qSg4sSuySahQxoAQAS664JbUvjCZLqHMzXMm3x7VPD5QZNWzOUeFVqQiGZkyrtPAACAPtaf1c6t5NemsPb8dT9n6tOh9F8q+DXJjIKh8s8Iy04KM7bxI3yp7BB6fCBthx63ToYWIQCAWe2/on1zL19wy+iYod7UgiFMZinvLCPZIfTaZPbEVFlWEJ17VAj4UZO9mV97RriOmRWdBi1CAACzesSNWjeadbp9zaJh3lTR6JvTJ7rYkewQOjuErleTwhrh4yrhD7/w47rSWT2o0QG0DE0Yk8LpBAAwK5oiTvdauc/lXvcj8FCQ7BC6ZAx7YCIb7UHN2y/4/nCjjahBG9FEEIQAABLg70DlRNLlaeze3xPRb8mNRMR6NEZCEAIASEnX3xNxT/rvbcQlmuzNfGmtVtD3rnpwGwQhAIAkBTreSMRdE9hoD2pOJd/tJy6ngi+/oPdtZoEQgsEyAABS192JyomkciLpk03aNae1f93OX1SRjO5UVhA9zAf3xHo4BCEAgJUI+j0Rqxq0y04Kz5XzLRxJ74ZEfAgEIQCAtYlwoyLcmDlRRJeIz5Txap5M6EZlBdGxPgjEOyEIAQCs1h2J+FQZzwkkLZDKDqGjPZCINyAIAQCs3x2JOHUjL6dJVg9qSg86zNXWExFBCABgQ+5IxNT1vIIhWT2oacF0LxcbTUQEIQCALdIl4hv9ybY67bKTwqgC3lVOsnpQjwXTITaWiAhCAADbRVMk1oeK9WE+irmRiCPyOXc7KqsH9XhPuqezTSQighAAAG5NREaXiLFruS52VFYPKjuE7uFkzYmIIAQAgJvaE/HDGKaiTrvspDBsDRfsTGUF0Vk9KD+lFSYighAAAO6BuSsRo1bxPZ2prCB6SjDtYy92faaDtUYBAOBBdIm4YAhzeqrs1b70nnpt+HJN7FpuwSGhTiV2caaAFiEAAOjFjiFpgXRaIFHzTEmtsOyEdk6lJsKNygqipwXTXpJtIyIIAQDAMIq7EvFflZpINyoriH6sJ+2hELs+AyEIAQCgg9oTUcUxpeeEZSe0b+7VDPWmsoLo9O60s0zs+vSDIAQAAGPZs3cm4sztNxIxozvtZNmJiCAEAACTaU/ExuvM6tPCspPCjG38CF8qK4ie2J12tMhERBACAIDpuchJdgidHUJfvU7W3J6Ik4JoB0sKH0uqBQAArI7r74nY0EbWnrktETODaKUFpJAFlAAAADbAze5GIl5pI/lnhGUnhb9t58d2pbN6UCkBtFy8ae0IQgAAMCv33xPxchspOCMsOiw8tYUfEyBaImJlGQAAEEcXO5IdQq9NZvdlsNEe1Lz9gt8PmuzN/Nozgka4+bLfGrVv7RXueK+aJ3/exrfxJigDQQgAACLzd6ByIunyNHbfxBuJ6Pt7InICCXWhLqi0r+y8GXrXBTJ5A/+IG2XHmODoCEIAALAUAb8nYuXvbUTfJZrpm/mUAKpJQ3RZeF0gmaX82K7UH8NME2G4RggAABYn0JHKiaRyIunj17Q/n9S+sVuobyNdFOTUNVrFa8cF0qZKQYIWIQAAWLJgZ+offel9E9kNY5kJ3eiiWpomWhOmIEEQAgCAJAQ5Ufvqtf/uL/g7UrdeLzQeghAAACxd+3XBZ0OEj2OoFo6YMAsRhAAAYNHuGB1DEfLJUMaEWYggBAAAi1bbop3S47Yxoros9LanbH0eYW1t7fXr18WuQkrOnTvX1tYmdhVScuHCBZVKJXYVUlJXV9fa2ip2FVJy6dKl5uZmsauwdEFO1GM9b6RVfX19U1MTIYQi5O+P0LY+jzArK+vw4cNiVyElTzzxxN69e8WuQkqeffbZbdu2iV2FlMyYMWPDhg1iVyElL7/8ckFBgdhVSMns2bOXL19u2n1KOAgBAACMhyAEAACbhiAEAACbRmm1WrFruOHy5cupqakXL17U8/UNDQ1OTk4si1Xi9IUzZqirV686ODjIZDKxC5GMxsZGe3t7uVwudiGSce3aNTs7Ozs7O7ELkYympiaZTKZQKPR8fX5+flhY2INfY0FBSAi5dOmSbjgQAACA8QICAh76y8yyghAAAMDMcI0QAABsGoIQAABsGoIQAABsGoIQAABsmiRH0tfV1e3evbu2tjYhISE4OFjsciSgurp63bp1tbW1QUFB2dnZzs7OYldk6TZt2lReXt7Q0BAQEPDEE094enqKXZFkLF68uEuXLmPGjBG7EEtXXl5eXV3d/vS5554TsRipuHr16nfffXfq1Cnd/5heXl4m2a0kW4QjRox45513Xn311V27doldizSkpKScOHEiMDBw/fr1/fv3b2xsFLsiS/fzzz9zHNejR49du3b16dOnrq5O7Iqk4fvvv//zn//88ccfi12IBCxZsuSHH3448Tuxy5GA8+fP9+/fv7y8vHv37rW1tUeOHDHVniU5fUIQBJqm+/XrN2vWrKlTp4pdjgSo1Wrd/FNBEEJCQt59993JkyeLXZRkhIeHz549e9q0aWIXYukuXLiQkJAwceLE3bt3FxUViV2OpXvhhRf8/f1fe+01sQuRjOnTp8tksq+++srke5Zki5CmJVm2iG5dhUGtVjs6OopYjLQcP3784sWLERERYhciAS+++OLcuXPd3NzELkQy9uzZ8/777y9duhT3R9NHYWFhVlZWbm7up59+euzYMRPuGYliW9555x1vb+/k5GSxC5GAN954IyAgIDw8/M033+zTp4/Y5Vi6n3/+Wa1WZ2Zmil2IZAQEBPj4+Fy9enX+/PlRUVHXrl0TuyKL1tTUVF9f//LLL1dVVZ04cWLgwIFbt2412d61ktW3b98ff/xR7CqkJDc319fX9/Dhw2IXIg0tLS0XLlxYuXJlly5dysrKxC7HotXX14eEhNTU1Gi12vnz548ePVrsiqSE47gBAwa8//77Yhdi0XS3L37nnXd0T19//fVx48aZaueSHDUKHbB06dJZs2Zt2LChV69eYtciDUqlUqlUZmRk5Ofnr1y5MjY2VuyKLFdxcXF9fX16ejohpK6u7tq1ayNGjNiyZYvYdUkDwzAxMTEYL/NgDg4Obm5u4eHhuqcRERF5eXmm2jm6Rm3CypUrZ86cuX79+oeuwg6EEJ7nr1+/rnus0Wj27dsXGBgobkkWLiUlpaSkZOHChQsXLszKyoqIiMDA0YdSqVS6B62traWlpbgO/VDp6enbt2/XPa6oqDDhGZPkqNEZM2ZUVFRUV1f7+Pi4u7t/8cUXAwYMELsoy6XRaBwcHDw8PPz8/HRbZsyYMX36dHGrsmT19fXh4eFDhw51cnLaunVrQEDAunXrlEql2HVJw4cfflhSUoJRow/l7+8fFRXl4uLyyy+/9O7dOz8/HzdjerDjx4+PHDkyNjaW47gdO3Zs2LAhNDTUJHuWZBAePXr01gvLoaGhTk5OItZj4bRabWVl5a1b/P39fXx8xKpHEs6ePVtZWalSqXr27BkdHS12OVKi6xoNCQkRuxBLV1NTU1lZqVar8RnT37Vr1zZs2KBQKIYNG2bChUEkGYQAAACmgmuEAABg0xCEAABg0xCEAABg0xCEAABg0xCEAABg0xCEAABg0xCEAABg0xCEAABg0xCEAABg0xCEAABg0xCEAABg0/4fQlHxlUSOFVIAAAAASUVORK5CYII=", "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", "\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", "\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.3" }, "kernelspec": { "name": "julia-1.6", "display_name": "Julia 1.6.3", "language": "julia" } }, "nbformat": 4 }