Elastic constants

We compute clamped-ion elastic constants of a crystal using the algorithmic differentiation density-functional perturbation theory (AD-DFPT) approach as introduced in [SPH25].

We consider a crystal in its equilibrium configuration, where all atomic forces and stresses vanish. Homogeneous strains $η$ are then applied relative to this relaxed structure. The elastic constants are derived from the stress-strain relationship. In Voigt notation, the stress $\sigma$ and strain $\eta$ tensors are represented as 6-component vectors. The elastic constants $C$ are then given by the Jacobian of the stress with respect to strain, forming a $6 \times 6$ matrix

\[ C = \frac{\partial \sigma}{\partial \eta}.\]

The sparsity pattern of the matrix $C$ follows from crystal symmetry and is tabulated in standard references (eg. Table 9 in [Nye1985]). This sparsity can be used a priori to reduce the number of strain patterns that need to be probed to extract all independent components of $C$. For example, cubic crystals have only three independent elastic constants $C_{11}$, $C_{12}$ and $C_{44}$, with the pattern

\[C = \begin{pmatrix} C_{11} & C_{12} & C_{12} & 0 & 0 & 0 \\ C_{12} & C_{11} & C_{12} & 0 & 0 & 0 \\ C_{12} & C_{12} & C_{11} & 0 & 0 & 0 \\ 0 & 0 & 0 & C_{44} & 0 & 0 \\ 0 & 0 & 0 & 0 & C_{44} & 0 \\ 0 & 0 & 0 & 0 & 0 & C_{44} \\ \end{pmatrix}.\]

Thus we can just choose a suitable strain pattern $\dot{\eta} = (1, 0, 0, 1, 0, 0)^\top$, such that $C\dot{\eta} = (C_{11}, C_{12}, C_{12}, C_{44}, 0, 0)^\top$. That is, for cubic crystals like diamond silicon we obtain all independent elastic constants from a single Jacobian-vector product on the stress-strain function.

This example computes the clamped-ion elastic tensor, keeping internal atomic positions fixed under strain. The relaxed-ion tensor includes additional corrections from internal relaxations, which can be obtained from first-order atomic displacements in DFPT (see [Wu2005]).

using DFTK
using PseudoPotentialData
using LinearAlgebra
using ForwardDiff
using DifferentiationInterface
using AtomsBuilder
using Unitful
using UnitfulAtomic

Computing PBE elastic constants

We start with the PBE [PBE1996] functional.

pseudopotentials = PseudoFamily("dojo.nc.sr.pbe.v0_4_1.standard.upf")
a0_pbe = 10.33u"bohr"  # Equilibrium lattice constant of silicon with PBE
model0 = model_DFT(bulk(:Si; a=a0_pbe); pseudopotentials, functionals=PBE())

Ecut = recommended_cutoff(model0).Ecut
kgrid = [4, 4, 4]
tol = 1e-6

function symmetries_from_strain(model0, voigt_strain)
    lattice = DFTK.voigt_strain_to_full(voigt_strain) * model0.lattice
    model = Model(model0; lattice, symmetries=true)
    model.symmetries
end

strain_pattern = [1., 0., 0., 1., 0., 0.];  # recovers [c11, c12, c12, c44, 0, 0]

For elastic constants beyond the bulk modulus, symmetry-breaking strains are required. That is, the symmetry group of the crystal is reduced. Here we simply precompute the relevant subgroup by applying the automatic symmetry detection (spglib) to the finitely perturbed crystal.

symmetries_strain = symmetries_from_strain(model0, 0.01 * strain_pattern)


function stress_from_strain(model0, voigt_strain; symmetries, Ecut, kgrid, tol)
    lattice = DFTK.voigt_strain_to_full(voigt_strain) * model0.lattice
    model = Model(model0; lattice, symmetries)
    basis = PlaneWaveBasis(model; Ecut, kgrid)
    scfres = self_consistent_field(basis; tol)
    DFTK.full_stress_to_voigt(compute_stresses_cart(scfres))
end

stress_fn(voigt_strain) = stress_from_strain(model0, voigt_strain;
                                             symmetries=symmetries_strain,
                                             Ecut, kgrid, tol)
stress, (dstress,) = value_and_pushforward(stress_fn, AutoForwardDiff(),
                                           zeros(6), (strain_pattern,));
n     Energy            log10(ΔE)   log10(Δρ)   Diag   Δtime 
---   ---------------   ---------   ---------   ----   ------
  1   -8.453495592272                   -0.94    5.3    416ms
  2   -8.455630031252       -2.67       -1.77    1.0    168ms
  3   -8.455778742205       -3.83       -2.89    1.9    198ms
  4   -8.455790083248       -4.95       -3.32    2.9    272ms
  5   -8.455790178610       -7.02       -3.95    1.2    162ms
  6   -8.455790187034       -8.07       -4.74    1.6    189ms
  7   -8.455790187698       -9.18       -5.67    2.1    221ms
  8   -8.455790187713      -10.82       -6.07    2.8    270ms
Solving response problem
Iter  Restart  Krydim  log10(res)  avg(CG)  Δtime   Comment
----  -------  ------  ----------  -------  ------  ---------------
                                      17.1   7.64s  Non-interacting
   1        0       1        0.11     12.7   15.7s  
   2        0       2       -0.70     11.4   826ms  
   3        0       3       -1.71     10.1   460ms  
   4        0       4       -2.78      8.1   387ms  
   5        0       5       -4.04      5.3   303ms  
   6        0       6       -5.57      1.9   190ms  
   7        0       7       -7.77      1.0   161ms  
   8        1       1       -6.38     16.3   923ms  Restart
                                      15.9   1.92s  Final orbitals

We can inspect the stress to verify it is small (close to equilibrium):

stress
6-element StaticArraysCore.SVector{6, Float64} with indices SOneTo(6):
 -2.3120131645660768e-5
 -2.312034974362303e-5
 -2.3120349743623055e-5
  7.034973824820316e-10
  0.0
  0.0

The response of the stress to strain_pattern contains the elastic constants in atomic units, with the expected pattern $(c11, c12, c12, c44, 0, 0)$:

dstress
6-element StaticArraysCore.SVector{6, Float64} with indices SOneTo(6):
 0.005319939387016797
 0.002024899420222948
 0.002024899420222948
 0.0033519068894872434
 0.0
 0.0

Convert to GPa:

println("C11: ", uconvert(u"GPa", dstress[1] * u"hartree" / u"bohr"^3))
println("C12: ", uconvert(u"GPa", dstress[2] * u"hartree" / u"bohr"^3))
println("C44: ", uconvert(u"GPa", dstress[4] * u"hartree" / u"bohr"^3))
C11: 156.5180202105481 GPa
C12: 59.5745976264789 GPa
C44: 98.61650520925065 GPa

These results can be compared directly to finite differences of the stress_fn:

h = 1e-3
dstress_fd = (stress_fn(h * strain_pattern) - stress_fn(-h * strain_pattern)) / 2h
println("C11 (FD): ", uconvert(u"GPa", dstress_fd[1] * u"hartree" / u"bohr"^3))
println("C12 (FD): ", uconvert(u"GPa", dstress_fd[2] * u"hartree" / u"bohr"^3))
println("C44 (FD): ", uconvert(u"GPa", dstress_fd[4] * u"hartree" / u"bohr"^3))
n     Energy            log10(ΔE)   log10(Δρ)   Diag   Δtime 
---   ---------------   ---------   ---------   ----   ------
  1   -8.453517556096                   -0.94    5.4    398ms
  2   -8.455633436047       -2.67       -1.77    1.0    156ms
  3   -8.455782741165       -3.83       -2.89    2.0    204ms
  4   -8.455795209965       -4.90       -3.31    2.9    282ms
  5   -8.455795304904       -7.02       -3.84    1.3    175ms
  6   -8.455795314288       -8.03       -4.92    1.4    177ms
  7   -8.455795315190       -9.04       -5.51    2.9    277ms
  8   -8.455795315197      -11.12       -5.94    1.7    202ms
  9   -8.455795315198      -12.21       -6.62    1.3    171ms
n     Energy            log10(ΔE)   log10(Δρ)   Diag   Δtime 
---   ---------------   ---------   ---------   ----   ------
  1   -8.453473298957                   -0.94    5.4    381ms
  2   -8.455620637872       -2.67       -1.78    1.0    154ms
  3   -8.455767643459       -3.83       -2.88    1.8    196ms
  4   -8.455782121422       -4.84       -3.25    2.9    274ms
  5   -8.455782238180       -6.93       -3.74    1.1    157ms
  6   -8.455782253233       -7.82       -4.79    1.4    181ms
  7   -8.455782254356       -8.95       -5.45    2.7    336ms
  8   -8.455782254370      -10.86       -5.94    1.8    847ms
  9   -8.455782254371      -12.04       -6.81    1.8    221ms
C11 (FD): 156.2909842093387 GPa
C12 (FD): 59.51630198935152 GPa
C44 (FD): 98.54130855707672 GPa

Here are AD-DFPT results from increasing discretization parameters:

Ecutkgridc11c12c44
18[4, 4, 4]156.5159.5798.61
18[8, 8, 8]153.5356.90100.07
24[8, 8, 8]153.2656.8299.97
24[14, 14, 14]153.0356.71100.09

For comparison, Materials Project for PBE relaxed-ion elastic constants of silicon mp-149: c11 = 153 GPa, c12 = 57 GPa, c44 = 74 GPa. Note the discrepancy in c44, which is due to us not yet including ionic relaxation in this example.

Moving to meta-GGA functionals

To make the problem a little more interesting we will now compute the elastic constants using the Laplacian-dependent r2SCAN-L functional (a deorbitalized version of r2SCAN, see [MT2020]), again using AD-DFPT.

model_r2scanl = model_DFT(bulk(:Si; a=a0_pbe);
                       pseudopotentials, functionals=[:mgga_x_r2scanl, :mgga_c_r2scanl])

stress, (dstress,) = value_and_pushforward(AutoForwardDiff(), zeros(6),
                                           (strain_pattern,)) do voigt_strain
    stress_from_strain(model_r2scanl, voigt_strain;
                       symmetries=symmetries_strain, Ecut, kgrid, tol)
end;
┌ Warning: Meta-GGAs with a Δρ term have not yet been thoroughly tested.
└ @ DFTK ~/work/DFTK.jl/DFTK.jl/src/terms/xc.jl:121
n     Energy            log10(ΔE)   log10(Δρ)   Diag   Δtime 
---   ---------------   ---------   ---------   ----   ------
  1   -8.490338765098                   -0.92    5.4    2.58s
  2   -8.493545534233       -2.49       -1.80    1.0    549ms
  3   -8.493817176106       -3.57       -2.63    2.2    264ms
  4   -8.493840375441       -4.63       -3.35    2.5    314ms
  5   -8.493846293412       -5.23       -3.74    2.3    282ms
  6   -8.493848958024       -5.57       -4.03    2.1    266ms
  7   -8.493850359820       -5.85       -4.37    2.8    310ms
  8   -8.493850688130       -6.48       -4.76    2.8    309ms
  9   -8.493850749962       -7.21       -5.18    3.1    327ms
 10   -8.493850753653       -8.43       -5.60    2.1    270ms
 11   -8.493850755126       -8.83       -5.91    2.2    284ms
 12   -8.493850755533       -9.39       -6.23    2.1    275ms
Solving response problem
Iter  Restart  Krydim  log10(res)  avg(CG)  Δtime   Comment
----  -------  ------  ----------  -------  ------  ---------------
                                      17.2   631ms  Non-interacting
   1        0       1        0.11     12.6   2.57s  
   2        0       2       -0.67     10.6   532ms  
   3        0       3       -1.30      8.4   459ms  
   4        0       4       -1.86      6.6   394ms  
   5        0       5       -2.56      5.0   345ms  
   6        0       6       -2.89      2.4   264ms  
   7        0       7       -3.35      2.2   257ms  
   8        0       8       -3.79      2.0   257ms  
   9        0       9       -4.26      1.9   249ms  
  10        0      10       -4.87      1.5   231ms  
  11        0      11       -7.04      1.0   246ms  
  12        1       1       -4.61     18.7   1.75s  Restart
  13        1       2       -4.97      4.6   318ms  
  14        1       3       -5.23      2.9   330ms  
  15        1       4       -5.60      3.1   288ms  
  16        1       5       -6.03      2.0   250ms  
  17        1       6       -6.43      2.0   251ms  
                                      16.5   713ms  Final orbitals

Note here that the value_and_pushforward(...) do voigt_strain ... end syntax defines an anonymous function to compute the stress from a given strain, i.e. defining the same kind of function as stress_fn in the above example. See the Julia documentation on the do-block-syntax for details.

We again inspect the r2SCAN-L stress, which is still small (despite using the PBE lattice constant), so we remain close to the equilibrium:

stress
6-element StaticArraysCore.SVector{6, Float64} with indices SOneTo(6):
 -1.991259797978997e-5
 -1.9912619197159416e-5
 -1.991261919715943e-5
  1.0278827069688885e-9
  0.0
  0.0

Finally, the r2SCAN-L elastic constants, which agree well with PBE:

println("C11: ", uconvert(u"GPa", dstress[1] * u"hartree" / u"bohr"^3))
println("C12: ", uconvert(u"GPa", dstress[2] * u"hartree" / u"bohr"^3))
println("C44: ", uconvert(u"GPa", dstress[4] * u"hartree" / u"bohr"^3))
C11: 155.91235129444064 GPa
C12: 59.34778650762075 GPa
C44: 97.94754812053009 GPa