Saving SCF results on disk and SCF checkpoints
For longer DFT calculations it is pretty standard to run them on a cluster in advance and to perform postprocessing (band structure calculation, plotting of density, etc.) at a later point and potentially on a different machine.
To support such workflows DFTK offers the two functions save_scfres and load_scfres, which allow to save the data structure returned by self_consistent_field on disk or retrieve it back into memory, respectively. For this purpose DFTK uses the JLD2.jl file format and Julia package. For the moment this process is considered an experimental feature and has a number of caveats, see the warnings below.
The load_scfres and save_scfres pair of functions are experimental features. This means:
- The interface of these functions as well as the format in which the data is stored on disk can change incompatibly in the future. At this point we make no promises ...
- JLD2 is not yet completely matured and it is recommended to only use it for short-term storage and not to archive scientific results.
- If you are using the functions to transfer data between different machines ensure that you use the same version of Julia, JLD2 and DFTK for saving and loading data.
To illustrate the use of the functions in practice we will compute the total energy of the O₂ molecule at PBE level. To get the triplet ground state we use a collinear spin polarisation (see Collinear spin and magnetic systems for details) and a bit of temperature to ease convergence:
using DFTK
using LinearAlgebra
using JLD2
d = 2.079  # oxygen-oxygen bondlength
a = 9.0    # size of the simulation box
lattice = diagm(a * ones(3))
O = ElementPsp(:O, psp=load_psp("hgh/pbe/O-q6.hgh"))
atoms = [O => d / 2a * [[0, 0, 1], [0, 0, -1]]]
magnetic_moments = [O => [1., 1.]]
Ecut  = 10  # Far too small to be converged
model = model_PBE(lattice, atoms, temperature=0.02, smearing=smearing=Smearing.Gaussian(),
                  magnetic_moments=magnetic_moments)
basis = PlaneWaveBasis(model; Ecut, kgrid=[1, 1, 1])
scfres = self_consistent_field(basis, tol=1e-2, ρ=guess_density(basis, magnetic_moments))
save_scfres("scfres.jld2", scfres);n Free energy Eₙ-Eₙ₋₁ ρout-ρin Magnet α Diag --- --------------- --------- -------- ------ ---- ---- 1 -27.63747131447 NaN 7.35e-01 0.001 0.80 5.0 2 -28.50135558247 -8.64e-01 4.82e-01 0.908 0.80 5.0 3 -28.90877560524 -4.07e-01 1.30e-01 1.279 0.80 5.0 4 -28.93724400390 -2.85e-02 6.85e-02 1.738 0.80 3.0 5 -28.93836789187 -1.12e-03 4.54e-02 1.933 0.80 3.0
scfres.energiesEnergy breakdown:
    Kinetic             16.9260762
    AtomicLocal         -58.8326844
    AtomicNonlocal      4.7447423 
    Ewald               -4.8994689
    PspCorrection       0.0044178 
    Hartree             19.5437626
    Xc                  -6.4228466
    Entropy             -0.0023669
    total               -28.938367891867
The scfres.jld2 file could now be transfered to a different computer, Where one could fire up a REPL to inspect the results of the above calculation:
using DFTK
using JLD2
loaded = load_scfres("scfres.jld2")
propertynames(loaded)(:ham, :basis, :energies, :converged, :ρ, :eigenvalues, :occupation, :εF, :n_iter, :n_ep_extra, :ψ, :diagonalization, :stage, :algorithm)
loaded.energiesEnergy breakdown:
    Kinetic             16.9260762
    AtomicLocal         -58.8326844
    AtomicNonlocal      4.7447423 
    Ewald               -4.8994689
    PspCorrection       0.0044178 
    Hartree             19.5437626
    Xc                  -6.4228466
    Entropy             -0.0023669
    total               -28.938367891867
Since the loaded data contains exactly the same data as the scfres returned by the SCF calculation one could use it to plot a band structure, e.g. plot_bandstructure(load_scfres("scfres.jld2")) directly from the stored data.
Checkpointing of SCF calculations
A related feature, which is very useful especially for longer calculations with DFTK is automatic checkpointing, where the state of the SCF is periodically written to disk. The advantage is that in case the calculation errors or gets aborted due to overrunning the walltime limit one does not need to start from scratch, but can continue the calculation from the last checkpoint.
To enable automatic checkpointing in DFTK one needs to pass the ScfSaveCheckpoints callback to self_consistent_field, for example:
callback = DFTK.ScfSaveCheckpoints()
scfres = self_consistent_field(basis; ρ=guess_density(basis, magnetic_moments),
                               tol=1e-2, callback=callback);Notice that using this callback makes the SCF go silent since the passed callback parameter overwrites the default value (namely DefaultScfCallback()) which exactly gives the familiar printing of the SCF convergence. If you want to have both (printing and checkpointing) you need to chain both callbacks:
callback = DFTK.ScfDefaultCallback() ∘ DFTK.ScfSaveCheckpoints(keep=true)
scfres = self_consistent_field(basis; ρ=guess_density(basis, magnetic_moments),
                               tol=1e-2, callback=callback);n Free energy Eₙ-Eₙ₋₁ ρout-ρin Magnet α Diag --- --------------- --------- -------- ------ ---- ---- 1 -27.63016482928 NaN 7.37e-01 0.001 0.80 5.5 2 -28.49928627351 -8.69e-01 4.83e-01 0.896 0.80 5.0 3 -28.90841916006 -4.09e-01 1.30e-01 1.266 0.80 5.0 4 -28.93715079108 -2.87e-02 6.89e-02 1.730 0.80 3.0 5 -28.93834951218 -1.20e-03 4.63e-02 1.930 0.80 3.0
For more details on using callbacks with DFTK's self_consistent_field function see Monitoring self-consistent field calculations.
By default checkpoint is saved in the file dftk_scf_checkpoint.jld2, which is deleted automatically once the SCF completes successfully. If one wants to keep the file one needs to specify keep=true as has been done in the ultimate SCF for demonstration purposes: now we can continue the previous calculation from the last checkpoint as if the SCF had been aborted. For this one just loads the checkpoint with load_scfres:
oldstate = load_scfres("dftk_scf_checkpoint.jld2")
scfres   = self_consistent_field(oldstate.basis, ρ=oldstate.ρ,
                                 ψ=oldstate.ψ, tol=1e-3);n Free energy Eₙ-Eₙ₋₁ ρout-ρin Magnet α Diag --- --------------- --------- -------- ------ ---- ---- 1 -28.93869756827 NaN 2.10e-02 1.985 0.80 1.0 2 -28.93937502126 -6.77e-04 1.14e-02 1.982 0.80 1.0
As JLD2 is an optional dependency of DFTK these three functions are only available once one has both imported DFTK and JLD2 (using DFTK and using JLD2).
(Cleanup files generated by this notebook)
rm("dftk_scf_checkpoint.jld2")
rm("scfres.jld2")