{ "cells": [ { "cell_type": "markdown", "source": [ "# Saving SCF results on disk and SCF checkpoints\n", "\n", "For longer DFT calculations it is pretty standard to run them on a cluster\n", "in advance and to perform postprocessing (band structure calculation,\n", "plotting of density, etc.) at a later point and potentially on a different\n", "machine.\n", "\n", "To support such workflows DFTK offers the two functions `save_scfres`\n", "and `load_scfres`, which allow to save the data structure returned\n", "by `self_consistent_field` on disk or retrieve it back into memory,\n", "respectively. For this purpose DFTK uses the\n", "[JLD2.jl](https://github.com/JuliaIO/JLD2.jl) file format and Julia package.\n", "For the moment this process is considered an experimental feature and\n", "has a number of caveats, see the warnings below.\n", "\n", "!!! warning \"Saving `scfres` is experimental\"\n", " The `load_scfres` and `save_scfres` pair of functions\n", " are experimental features. This means:\n", "\n", " - The interface of these functions\n", " as well as the format in which the data is stored on disk can\n", " change incompatibly in the future. At this point we make no promises ...\n", " - JLD2 is not yet completely matured\n", " and it is recommended to only use it for short-term storage\n", " and **not** to archive scientific results.\n", " - If you are using the functions to transfer data between different\n", " machines ensure that you use the **same version of Julia, JLD2 and DFTK**\n", " for saving and loading data.\n", "\n", "To illustrate the use of the functions in practice we will compute\n", "the total energy of the O₂ molecule at PBE level. To get the triplet\n", "ground state we use a collinear spin polarisation\n", "(see Collinear spin and magnetic systems for details)\n", "and a bit of temperature to ease convergence:" ], "metadata": {} }, { "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "n Free energy Eₙ-Eₙ₋₁ ρout-ρin Magnet Diag\n", "--- --------------- --------- -------- ------ ----\n", " 1 -27.63593686918 NaN 9.78e-01 0.001 5.0 \n", " 2 -28.50059060049 -8.65e-01 6.72e-01 0.910 5.0 \n", " 3 -28.90870777240 -4.08e-01 1.52e-01 1.282 5.0 \n", " 4 -28.93726448188 -2.86e-02 3.94e-02 1.739 3.0 \n", " 5 -28.93850097313 -1.24e-03 2.90e-02 1.933 2.0 \n" ] } ], "cell_type": "code", "source": [ "using DFTK\n", "using LinearAlgebra\n", "using JLD2\n", "\n", "d = 2.079 # oxygen-oxygen bondlength\n", "a = 9.0 # size of the simulation box\n", "lattice = diagm(a * ones(3))\n", "O = ElementPsp(:O, psp=load_psp(\"hgh/pbe/O-q6.hgh\"))\n", "atoms = [O => d / 2a * [[0, 0, 1], [0, 0, -1]]]\n", "magnetic_moments = [O => [1., 1.]]\n", "\n", "Ecut = 10 # Far too small to be converged\n", "model = model_PBE(lattice, atoms, temperature=0.02, smearing=smearing=Smearing.Gaussian(),\n", " magnetic_moments=magnetic_moments)\n", "basis = PlaneWaveBasis(model, Ecut; kgrid=[1, 1, 1])\n", "\n", "ρspin = guess_spin_density(basis, magnetic_moments)\n", "scfres = self_consistent_field(basis, tol=1e-2, ρspin=ρspin)\n", "save_scfres(\"scfres.jld2\", scfres);" ], "metadata": {}, "execution_count": 1 }, { "outputs": [ { "output_type": "execute_result", "data": { "text/plain": "Energy breakdown:\n Kinetic 16.9137119\n AtomicLocal -58.8113237\n AtomicNonlocal 4.7461213 \n Ewald -4.8994689\n PspCorrection 0.0044178 \n Hartree 19.5308889\n Xc -6.4204952\n Entropy -0.0023531\n\n total -28.938500973131\n" }, "metadata": {}, "execution_count": 2 } ], "cell_type": "code", "source": [ "scfres.energies" ], "metadata": {}, "execution_count": 2 }, { "cell_type": "markdown", "source": [ "The `scfres.jld2` file could now be transfered to a different computer,\n", "Where one could fire up a REPL to inspect the results of the above\n", "calculation:" ], "metadata": {} }, { "outputs": [ { "output_type": "execute_result", "data": { "text/plain": "(:ham, :basis, :energies, :converged, :ρ, :ρspin, :eigenvalues, :occupation, :εF, :n_iter, :n_ep_extra, :ψ, :diagonalization, :stage)" }, "metadata": {}, "execution_count": 3 } ], "cell_type": "code", "source": [ "using DFTK\n", "using JLD2\n", "loaded = load_scfres(\"scfres.jld2\")\n", "propertynames(loaded)" ], "metadata": {}, "execution_count": 3 }, { "outputs": [ { "output_type": "execute_result", "data": { "text/plain": "Energy breakdown:\n Kinetic 16.9137119\n AtomicLocal -58.8113237\n AtomicNonlocal 4.7461213 \n Ewald -4.8994689\n PspCorrection 0.0044178 \n Hartree 19.5308889\n Xc -6.4204952\n Entropy -0.0023531\n\n total -28.938500973131\n" }, "metadata": {}, "execution_count": 4 } ], "cell_type": "code", "source": [ "loaded.energies" ], "metadata": {}, "execution_count": 4 }, { "cell_type": "markdown", "source": [ "Since the loaded data contains exactly the same data as the `scfres` returned by the\n", "SCF calculation one could use it to plot a band structure, e.g.\n", "`plot_bandstructure(load_scfres(\"scfres.jld2\"))` directly from the stored data." ], "metadata": {} }, { "cell_type": "markdown", "source": [ "## Checkpointing of SCF calculations\n", "A related feature, which is very useful especially for longer calculations with DFTK\n", "is automatic checkpointing, where the state of the SCF is periodically written to disk.\n", "The advantage is that in case the calculation errors or gets aborted due\n", "to overrunning the walltime limit one does not need to start from scratch,\n", "but can continue the calculation from the last checkpoint.\n", "\n", "To enable automatic checkpointing in DFTK one needs to pass the `ScfSaveCheckpoints`\n", "callback to `self_consistent_field`, for example:" ], "metadata": {} }, { "outputs": [], "cell_type": "code", "source": [ "callback = DFTK.ScfSaveCheckpoints()\n", "scfres = self_consistent_field(basis, tol=1e-2, ρspin=ρspin, callback=callback);" ], "metadata": {}, "execution_count": 5 }, { "cell_type": "markdown", "source": [ "Notice that using this callback makes the SCF go silent since the passed\n", "callback parameter overwrites the default value (namely `DefaultScfCallback()`)\n", "which exactly gives the familiar printing of the SCF convergence.\n", "If you want to have both (printing and checkpointing) you need to chain\n", "both callbacks:" ], "metadata": {} }, { "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "n Free energy Eₙ-Eₙ₋₁ ρout-ρin Magnet Diag\n", "--- --------------- --------- -------- ------ ----\n", " 1 -27.63578746041 NaN 9.78e-01 0.001 5.0 \n", " 2 -28.50089669139 -8.65e-01 6.72e-01 0.904 5.0 \n", " 3 -28.90866802881 -4.08e-01 1.52e-01 1.275 5.0 \n", " 4 -28.93721740399 -2.85e-02 3.94e-02 1.735 3.0 \n", " 5 -28.93848884918 -1.27e-03 2.91e-02 1.932 2.0 \n" ] } ], "cell_type": "code", "source": [ "callback = DFTK.ScfDefaultCallback() ∘ DFTK.ScfSaveCheckpoints(keep=true)\n", "scfres = self_consistent_field(basis, tol=1e-2, ρspin=ρspin, callback=callback);" ], "metadata": {}, "execution_count": 6 }, { "cell_type": "markdown", "source": [ "For more details on using callbacks with DFTK's `self_consistent_field` function\n", "see Monitoring self-consistent field calculations." ], "metadata": {} }, { "cell_type": "markdown", "source": [ "By default checkpoint is saved in the file `dftk_scf_checkpoint.jld2`, which is\n", "deleted automatically once the SCF completes successfully. If one wants to keep\n", "the file one needs to specify `keep=true` as has been done in the ultimate SCF\n", "for demonstration purposes: now we can continue the previous calculation\n", "from the last checkpoint as if the SCF had been aborted.\n", "For this one just loads the checkpoint with `load_scfres`:" ], "metadata": {} }, { "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "n Free energy Eₙ-Eₙ₋₁ ρout-ρin Magnet Diag\n", "--- --------------- --------- -------- ------ ----\n", " 1 -28.93886097033 NaN 2.67e-02 1.985 1.0 \n", " 2 -28.93943540799 -5.74e-04 1.30e-02 1.982 1.0 \n" ] } ], "cell_type": "code", "source": [ "oldstate = load_scfres(\"dftk_scf_checkpoint.jld2\")\n", "scfres = self_consistent_field(oldstate.basis, ρ=oldstate.ρ, ρspin=oldstate.ρspin,\n", " ψ=oldstate.ψ, tol=1e-3);" ], "metadata": {}, "execution_count": 7 }, { "cell_type": "markdown", "source": [ "!!! note \"Availability of `load_scfres`, `save_scfres` and `ScfSaveCheckpoints`\"\n", " As JLD2 is an optional dependency of DFTK these three functions are only\n", " available once one has *both* imported DFTK and JLD2 (`using DFTK`\n", " and `using JLD2`)." ], "metadata": {} }, { "cell_type": "markdown", "source": [ "(Cleanup files generated by this notebook)" ], "metadata": {} }, { "outputs": [], "cell_type": "code", "source": [ "rm(\"dftk_scf_checkpoint.jld2\")\n", "rm(\"scfres.jld2\")" ], "metadata": {}, "execution_count": 8 } ], "nbformat_minor": 3, "metadata": { "language_info": { "file_extension": ".jl", "mimetype": "application/julia", "name": "julia", "version": "1.6.0" }, "kernelspec": { "name": "julia-1.6", "display_name": "Julia 1.6.0", "language": "julia" } }, "nbformat": 4 }