Custom solvers

In this example, we show how to define custom solvers. Our system will again be silicon, because we are not very imaginative

using DFTK, LinearAlgebra

a = 10.26
lattice = a / 2 * [[0 1 1.];
                   [1 0 1.];
                   [1 1 0.]]
Si = ElementPsp(:Si; psp=load_psp("hgh/lda/Si-q4"))
atoms = [Si, Si]
positions =  [ones(3)/8, -ones(3)/8]

# We take very (very) crude parameters
model = model_DFT(lattice, atoms, positions; functionals=LDA())
basis = PlaneWaveBasis(model; Ecut=5, kgrid=[1, 1, 1]);

We define our custom fix-point solver: simply a damped fixed-point

function my_fp_solver(f, x0, info0; maxiter)
    mixing_factor = .7
    x = x0
    info = info0
    for n = 1:maxiter
        fx, info = f(x, info)
        if info.converged || info.timedout
            break
        end
        x = x + mixing_factor * (fx - x)
    end
    (; fixpoint=x, info)
end;

Note that the fixpoint map f operates on an auxiliary variable info for state bookkeeping. Early termination criteria are flagged from inside the function f using boolean flags info.converged and info.timedout. For control over these criteria, see the is_converged and maxtime keyword arguments of self_consistent_field.

Our eigenvalue solver just forms the dense matrix and diagonalizes it explicitly (this only works for very small systems)

function my_eig_solver(A, X0; maxiter, tol, kwargs...)
    n = size(X0, 2)
    A = Array(A)
    E = eigen(A)
    λ = E.values[1:n]
    X = E.vectors[:, 1:n]
    (; λ, X, residual_norms=[], n_iter=0, converged=true, n_matvec=0)
end;

Finally we also define our custom mixing scheme. It will be a mixture of simple mixing (for the first 2 steps) and than default to Kerker mixing. In the mixing interface δF is $(ρ_\text{out} - ρ_\text{in})$, i.e. the difference in density between two subsequent SCF steps and the mix function returns $δρ$, which is added to $ρ_\text{in}$ to yield $ρ_\text{next}$, the density for the next SCF step.

struct MyMixing
    n_simple  # Number of iterations for simple mixing
end
MyMixing() = MyMixing(2)

function DFTK.mix_density(mixing::MyMixing, basis, δF; n_iter, kwargs...)
    if n_iter <= mixing.n_simple
        return δF  # Simple mixing -> Do not modify update at all
    else
        # Use the default KerkerMixing from DFTK
        DFTK.mix_density(KerkerMixing(), basis, δF; kwargs...)
    end
end

That's it! Now we just run the SCF with these solvers

scfres = self_consistent_field(basis;
                               tol=1e-4,
                               solver=my_fp_solver,
                               eigensolver=my_eig_solver,
                               mixing=MyMixing());
n     Energy            log10(ΔE)   log10(Δρ)   Diag   Δtime
---   ---------------   ---------   ---------   ----   ------
  1   -7.090073978406                   -0.39    0.0    292ms
  2   -7.226157343674       -0.87       -0.65    0.0    330ms
  3   -7.249782177445       -1.63       -1.17    0.0   41.5ms
  4   -7.250986294028       -2.92       -1.48    0.0   52.5ms
  5   -7.251258312426       -3.57       -1.78    0.0   52.7ms
  6   -7.251319808416       -4.21       -2.07    0.0    106ms
  7   -7.251334099764       -4.84       -2.35    0.0   40.1ms
  8   -7.251337569219       -5.46       -2.62    0.0   45.5ms
  9   -7.251338457710       -6.05       -2.88    0.0   51.0ms
 10   -7.251338698748       -6.62       -3.14    0.0   92.2ms
 11   -7.251338767946       -7.16       -3.39    0.0   40.9ms
 12   -7.251338788853       -7.68       -3.64    0.0   40.8ms
 13   -7.251338795449       -8.18       -3.88    0.0   51.9ms
 14   -7.251338797603       -8.67       -4.12    0.0   71.8ms

Note that the default convergence criterion is the difference in density. When this gets below tol, the fixed-point solver terminates. You can also customize this with the is_converged keyword argument to self_consistent_field, as shown below.

Customizing the convergence criterion

Here is an example of a defining a custom convergence criterion and specifying it using the is_converged callback keyword to self_consistent_field.

function my_convergence_criterion(info)
    tol = 1e-10
    length(info.history_Etot) < 2 && return false
    ΔE = (info.history_Etot[end-1] - info.history_Etot[end])
    ΔE < tol
end

scfres2 = self_consistent_field(basis;
                                solver=my_fp_solver,
                                is_converged=my_convergence_criterion,
                                eigensolver=my_eig_solver,
                                mixing=MyMixing());
n     Energy            log10(ΔE)   log10(Δρ)   Diag   Δtime
---   ---------------   ---------   ---------   ----   ------
  1   -7.090073978406                   -0.39    0.0    119ms
  2   -7.226157343674       -0.87       -0.65    0.0    116ms
  3   -7.249782177445       -1.63       -1.17    0.0   69.5ms
  4   -7.250986294028       -2.92       -1.48    0.0   40.3ms
  5   -7.251258312426       -3.57       -1.78    0.0   48.6ms
  6   -7.251319808416       -4.21       -2.07    0.0   83.8ms
  7   -7.251334099764       -4.84       -2.35    0.0   39.9ms
  8   -7.251337569219       -5.46       -2.62    0.0   42.5ms
  9   -7.251338457710       -6.05       -2.88    0.0   50.9ms
 10   -7.251338698748       -6.62       -3.14    0.0   66.2ms
 11   -7.251338767946       -7.16       -3.39    0.0   40.9ms
 12   -7.251338788853       -7.68       -3.64    0.0   63.6ms
 13   -7.251338795449       -8.18       -3.88    0.0   39.1ms
 14   -7.251338797603       -8.67       -4.12    0.0   39.4ms
 15   -7.251338798325       -9.14       -4.36    0.0   49.4ms
 16   -7.251338798572       -9.61       -4.59    0.0   39.5ms
 17   -7.251338798658      -10.07       -4.82    0.0   51.3ms