Eigenvalues of the dielectric matrix

We compute a few eigenvalues of the dielectric matrix ($q=0$, $ω=0$) iteratively.

using DFTK
using Plots
using KrylovKit
using Printf

# Calculation parameters
kgrid = [1, 1, 1]
Ecut = 5

# Silicon lattice
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]

# Compute the dielectric operator without symmetries
model  = model_LDA(lattice, atoms, positions, symmetries=false)
basis  = PlaneWaveBasis(model; Ecut, kgrid)
scfres = self_consistent_field(basis, tol=1e-8);
n     Energy            log10(ΔE)   log10(Δρ)   Diag   Δtime
---   ---------------   ---------   ---------   ----   ------
  1   -7.236059688741                   -0.50    8.0    118ms
  2   -7.250580639155       -1.84       -1.41    1.0   7.63ms
  3   -7.251252635045       -3.17       -2.25    1.0   6.86ms
  4   -7.251298963481       -4.33       -2.31    2.0   6.23ms
  5   -7.251334791947       -4.45       -2.76    1.0   5.13ms
  6   -7.251338584444       -5.42       -3.42    1.0   5.19ms
  7   -7.251338756391       -6.76       -3.80    1.0   5.26ms
  8   -7.251338791516       -7.45       -4.13    2.0   6.20ms
  9   -7.251338798324       -8.17       -4.63    2.0   6.45ms
 10   -7.251338798575       -9.60       -4.98    2.0   6.19ms
 11   -7.251338798688       -9.95       -5.65    1.0   5.70ms
 12   -7.251338798702      -10.84       -5.76    3.0   7.68ms
 13   -7.251338798704      -11.73       -6.18    1.0   5.74ms
 14   -7.251338798704      -12.28       -6.73    2.0   6.69ms
 15   -7.251338798705      -13.35       -7.31    1.0   5.71ms
 16   -7.251338798705      -14.45       -7.26    3.0   7.39ms
 17   -7.251338798705      -14.57       -7.69    1.0   5.77ms
 18   -7.251338798705   +    -Inf       -8.08    2.0   6.66ms

Applying $ε^† ≔ (1- χ_0 K)$

function eps_fun(δρ)
    δV = apply_kernel(basis, δρ; ρ=scfres.ρ)
    χ0δV = apply_χ0(scfres, δV)
    δρ - χ0δV
end;

… eagerly diagonalizes the subspace matrix at each iteration

eigsolve(eps_fun, randn(size(scfres.ρ)), 5, :LM; eager=true, verbosity=3);
[ Info: Arnoldi iteration step 1: normres = 0.032740044888630186
[ Info: Arnoldi iteration step 2: normres = 0.6457263164858278
[ Info: Arnoldi iteration step 3: normres = 0.8651609394525643
[ Info: Arnoldi iteration step 4: normres = 0.39463105163081874
[ Info: Arnoldi iteration step 5: normres = 0.2154823579081172
[ Info: Arnoldi schursolve in iter 1, krylovdim = 5: 0 values converged, normres = (2.54e-02, 4.00e-02, 1.48e-01, 1.31e-01, 7.12e-02)
[ Info: Arnoldi iteration step 6: normres = 0.52073582234337
[ Info: Arnoldi schursolve in iter 1, krylovdim = 6: 0 values converged, normres = (9.47e-03, 4.40e-02, 5.03e-01, 1.04e-01, 2.80e-02)
[ Info: Arnoldi iteration step 7: normres = 0.13097193483628333
[ Info: Arnoldi schursolve in iter 1, krylovdim = 7: 0 values converged, normres = (7.44e-04, 1.49e-02, 6.73e-02, 4.81e-02, 9.53e-02)
[ Info: Arnoldi iteration step 8: normres = 0.09055720364947988
[ Info: Arnoldi schursolve in iter 1, krylovdim = 8: 0 values converged, normres = (2.95e-05, 9.91e-04, 4.91e-03, 1.57e-02, 5.10e-02)
[ Info: Arnoldi iteration step 9: normres = 0.07294893988575087
[ Info: Arnoldi schursolve in iter 1, krylovdim = 9: 0 values converged, normres = (9.32e-07, 5.18e-05, 2.85e-04, 4.08e-03, 2.54e-02)
[ Info: Arnoldi iteration step 10: normres = 0.09181858206809758
[ Info: Arnoldi schursolve in iter 1, krylovdim = 10: 0 values converged, normres = (3.68e-08, 3.36e-06, 2.05e-05, 1.21e-03, 1.50e-02)
[ Info: Arnoldi iteration step 11: normres = 0.0719246920669708
[ Info: Arnoldi schursolve in iter 1, krylovdim = 11: 0 values converged, normres = (1.15e-09, 1.72e-07, 1.17e-06, 2.92e-04, 7.99e-03)
[ Info: Arnoldi iteration step 12: normres = 0.08430572209295788
[ Info: Arnoldi schursolve in iter 1, krylovdim = 12: 0 values converged, normres = (4.12e-11, 1.01e-08, 7.61e-08, 7.38e-05, 4.11e-03)
[ Info: Arnoldi iteration step 13: normres = 0.0417222003335874
[ Info: Arnoldi schursolve in iter 1, krylovdim = 13: 1 values converged, normres = (7.40e-13, 2.98e-10, 2.49e-09, 9.95e-06, 1.22e-03)
[ Info: Arnoldi iteration step 14: normres = 0.615688904531361
[ Info: Arnoldi schursolve in iter 1, krylovdim = 14: 1 values converged, normres = (2.16e-13, 1.53e-10, 1.45e-09, 5.53e-05, 5.47e-01)
[ Info: Arnoldi iteration step 15: normres = 0.08802684161421216
[ Info: Arnoldi schursolve in iter 1, krylovdim = 15: 1 values converged, normres = (1.76e-14, 3.16e-10, 7.79e-02, 6.46e-04, 2.90e-05)
[ Info: Arnoldi iteration step 16: normres = 0.7134492142121792
[ Info: Arnoldi schursolve in iter 1, krylovdim = 16: 1 values converged, normres = (7.62e-15, 4.26e-10, 7.64e-02, 4.90e-04, 7.03e-01)
[ Info: Arnoldi iteration step 17: normres = 0.03800825582611777
[ Info: Arnoldi schursolve in iter 1, krylovdim = 17: 1 values converged, normres = (1.96e-16, 3.92e-09, 1.14e-02, 2.67e-06, 2.36e-02)
[ Info: Arnoldi iteration step 18: normres = 0.02450730348218596
[ Info: Arnoldi schursolve in iter 1, krylovdim = 18: 1 values converged, normres = (1.98e-18, 8.17e-09, 1.84e-04, 4.00e-04, 1.18e-04)
[ Info: Arnoldi iteration step 19: normres = 0.0494134037702408
[ Info: Arnoldi schursolve in iter 1, krylovdim = 19: 1 values converged, normres = (4.06e-20, 5.69e-06, 2.09e-06, 1.76e-07, 1.51e-05)
[ Info: Arnoldi iteration step 20: normres = 0.22862432452198292
[ Info: Arnoldi schursolve in iter 1, krylovdim = 20: 1 values converged, normres = (4.20e-21, 4.39e-09, 1.06e-06, 1.72e-06, 2.41e-06)
[ Info: Arnoldi iteration step 21: normres = 0.024879894428271408
[ Info: Arnoldi schursolve in iter 1, krylovdim = 21: 1 values converged, normres = (4.78e-23, 2.05e-08, 3.99e-09, 6.26e-08, 2.03e-08)
[ Info: Arnoldi iteration step 22: normres = 0.015611631081200405
[ Info: Arnoldi schursolve in iter 1, krylovdim = 22: 1 values converged, normres = (3.06e-25, 1.36e-10, 1.64e-10, 4.99e-10, 5.43e-10)
[ Info: Arnoldi iteration step 23: normres = 0.06563496120987931
[ Info: Arnoldi schursolve in iter 1, krylovdim = 23: 1 values converged, normres = (8.23e-27, 5.70e-12, 7.12e-12, 2.39e-11, 2.52e-11)
[ Info: Arnoldi iteration step 24: normres = 0.2812544128866648
[ Info: Arnoldi schursolve in iter 1, krylovdim = 24: 1 values converged, normres = (2.09e-27, 8.01e-12, 9.99e-12, 1.10e-10, 1.18e-10)
[ Info: Arnoldi iteration step 25: normres = 0.05142710060754561
[ Info: Arnoldi schursolve in iter 1, krylovdim = 25: 3 values converged, normres = (4.69e-29, 3.88e-13, 4.85e-13, 1.64e-08, 2.77e-08)
[ Info: Arnoldi iteration step 26: normres = 0.06745806197730123
[ Info: Arnoldi schursolve in iter 1, krylovdim = 26: 3 values converged, normres = (1.43e-30, 2.00e-14, 2.50e-14, 2.74e-09, 4.02e-09)
[ Info: Arnoldi iteration step 27: normres = 0.024260899698532647
[ Info: Arnoldi schursolve in iter 1, krylovdim = 27: 3 values converged, normres = (1.43e-32, 3.19e-16, 3.99e-16, 7.22e-10, 7.23e-10)
[ Info: Arnoldi iteration step 28: normres = 0.10018013434274649
[ Info: Arnoldi schursolve in iter 1, krylovdim = 28: 3 values converged, normres = (6.16e-34, 2.26e-17, 2.83e-17, 4.09e-08, 1.95e-10)
[ Info: Arnoldi iteration step 29: normres = 0.018611704430434035
[ Info: Arnoldi schursolve in iter 1, krylovdim = 29: 3 values converged, normres = (4.86e-36, 2.91e-19, 3.63e-19, 1.33e-09, 1.13e-08)
[ Info: Arnoldi iteration step 30: normres = 0.24079609721844034
[ Info: Arnoldi schursolve in iter 1, krylovdim = 30: 3 values converged, normres = (5.14e-37, 5.10e-20, 6.37e-20, 1.82e-09, 1.42e-09)
[ Info: Arnoldi schursolve in iter 2, krylovdim = 19: 3 values converged, normres = (5.14e-37, 5.10e-20, 6.37e-20, 1.82e-09, 1.42e-09)
[ Info: Arnoldi iteration step 20: normres = 0.04625973655445402
[ Info: Arnoldi schursolve in iter 2, krylovdim = 20: 3 values converged, normres = (1.12e-38, 1.95e-21, 2.44e-21, 7.67e-11, 6.44e-11)
[ Info: Arnoldi iteration step 21: normres = 0.07587243281295612
[ Info: Arnoldi schursolve in iter 2, krylovdim = 21: 3 values converged, normres = (3.60e-40, 1.03e-22, 1.28e-22, 4.46e-12, 3.74e-12)
[ Info: Arnoldi iteration step 22: normres = 0.015865725934302136
┌ Info: Arnoldi eigsolve finished after 2 iterations:
*  6 eigenvalues converged
*  norm of residuals = (2.4073050800734227e-42, 1.1087520092633694e-24, 1.039247252607569e-24, 5.337150052321752e-14, 7.18484904407011e-15, 1.2707570950656154e-15)
*  number of operations = 33