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 UnitfulAtomicComputing PBE elastic constants
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.453458239089 -0.94 5.4 419ms
2 -8.455625373835 -2.66 -1.77 1.0 151ms
3 -8.455778064722 -3.82 -2.88 2.0 206ms
4 -8.455790090604 -4.92 -3.33 2.9 274ms
5 -8.455790177610 -7.06 -3.87 1.2 165ms
6 -8.455790187019 -8.03 -4.91 1.4 180ms
7 -8.455790187701 -9.17 -5.34 2.9 345ms
8 -8.455790187712 -10.93 -5.95 1.6 807ms
9 -8.455790187713 -12.06 -6.63 1.7 190ms
Solving response problem
Iter Restart Krydim log10(res) avg(CG) Δtime Comment
---- ------- ------ ---------- ------- ------ ---------------
17.9 7.86s Non-interacting
1 0 1 0.11 13.7 16.0s
2 0 2 -0.70 12.4 857ms
3 0 3 -1.71 10.8 488ms
4 0 4 -2.78 8.9 497ms
5 0 5 -4.04 6.6 955ms
6 0 6 -5.56 2.7 208ms
7 0 7 -7.80 1.0 155ms
8 1 1 -6.72 17.8 962ms Restart
16.9 1.90s Final orbitals
We can inspect the stress to verify it is small (close to equilibrium):
stress6-element StaticArraysCore.SVector{6, Float64} with indices SOneTo(6):
-2.3121181250480202e-5
-2.312120174419833e-5
-2.3121201744198344e-5
-3.125323634164981e-10
0.0
0.0The response of the stress to strain_pattern contains the elastic constants in atomic units, with the expected pattern $(c11, c12, c12, c44, 0, 0)$:
dstress6-element StaticArraysCore.SVector{6, Float64} with indices SOneTo(6):
0.005319943636130809
0.0020248998197726267
0.0020248998197726267
0.003351909495927732
0.0
0.0Convert 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.51814522379811 GPa
C12: 59.57460938163626 GPa
C44: 98.61658189337719 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.453510463056 -0.94 5.3 418ms
2 -8.455640352514 -2.67 -1.77 1.0 162ms
3 -8.455784844910 -3.84 -2.88 1.9 213ms
4 -8.455795213411 -4.98 -3.38 3.0 990ms
5 -8.455795306845 -7.03 -3.99 1.4 169ms
6 -8.455795314861 -8.10 -5.22 1.8 195ms
7 -8.455795315197 -9.47 -5.84 3.1 278ms
8 -8.455795315198 -11.82 -6.67 1.6 177ms
n Energy log10(ΔE) log10(Δρ) Diag Δtime
--- --------------- --------- --------- ---- ------
1 -8.453471721038 -0.94 5.4 380ms
2 -8.455621978386 -2.67 -1.77 1.0 189ms
3 -8.455767558760 -3.84 -2.88 1.9 204ms
4 -8.455782124323 -4.84 -3.24 3.0 290ms
5 -8.455782237569 -6.95 -3.71 1.1 162ms
6 -8.455782252866 -7.82 -4.63 1.3 167ms
7 -8.455782254339 -8.83 -5.31 2.3 236ms
8 -8.455782254370 -10.50 -5.89 2.2 229ms
9 -8.455782254371 -11.91 -6.50 1.9 210ms
C11 (FD): 156.28380911427848 GPa
C12 (FD): 59.5081046827513 GPa
C44 (FD): 98.53544755238303 GPa
Here are AD-DFPT results from increasing discretization parameters:
| Ecut | kgrid | c11 | c12 | c44 |
|---|---|---|---|---|
| 18 | [4, 4, 4] | 156.51 | 59.57 | 98.61 |
| 18 | [8, 8, 8] | 153.53 | 56.90 | 100.07 |
| 24 | [8, 8, 8] | 153.26 | 56.82 | 99.97 |
| 24 | [14, 14, 14] | 153.03 | 56.71 | 100.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 constant using deorbitalised r2-SCAN functional, again using AD-DFPT.
model_scan = model_DFT(bulk(:Si; a=a0_pbe);
pseudopotentials, functionals=[:mgga_c_r2scanl, :mgga_x_r2scanl])
stress, (dstress,) = value_and_pushforward(AutoForwardDiff(), zeros(6),
(strain_pattern,)) do voigt_strain
stress_from_strain(model_scan, 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.490329113249 -0.92 5.6 2.67s
2 -8.493542336640 -2.49 -1.80 1.0 534ms
3 -8.493818689110 -3.56 -2.62 2.1 281ms
4 -8.493840490786 -4.66 -3.34 2.3 356ms
5 -8.493846265455 -5.24 -3.73 2.5 934ms
6 -8.493849055246 -5.55 -4.04 2.1 262ms
7 -8.493850391470 -5.87 -4.38 2.8 297ms
8 -8.493850696848 -6.52 -4.77 2.6 295ms
9 -8.493850750205 -7.27 -5.15 3.2 328ms
10 -8.493850753696 -8.46 -5.62 2.0 255ms
11 -8.493850755152 -8.84 -5.92 2.3 287ms
12 -8.493850755536 -9.42 -6.22 2.0 267ms
Solving response problem
Iter Restart Krydim log10(res) avg(CG) Δtime Comment
---- ------- ------ ---------- ------- ------ ---------------
17.2 623ms Non-interacting
1 0 1 0.11 12.6 2.56s
2 0 2 -0.67 10.6 1.27s
3 0 3 -1.30 8.4 443ms
4 0 4 -1.86 6.6 396ms
5 0 5 -2.56 5.0 334ms
6 0 6 -2.89 2.3 257ms
7 0 7 -3.35 2.2 256ms
8 0 8 -3.79 2.0 249ms
9 0 9 -4.26 1.9 254ms
10 0 10 -4.92 1.5 231ms
11 0 11 -7.09 1.0 223ms
12 1 1 -4.60 18.7 979ms Restart
13 1 2 -4.95 4.6 326ms
14 1 3 -5.22 2.8 280ms
15 1 4 -5.58 3.1 280ms
16 1 5 -6.02 2.0 257ms
17 1 6 -6.42 2.0 253ms
16.4 1.43s 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 stress, which is still small (despite using the PBE lattice constant), so we remain close to the equilibrium:
stress6-element StaticArraysCore.SVector{6, Float64} with indices SOneTo(6):
-1.991282002855195e-5
-1.9912843842193436e-5
-1.991284384219345e-5
4.036723586113369e-10
0.0
0.0Finially the r2-SCAN 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.91237678127987 GPa
C12: 59.347787423094516 GPa
C44: 97.94755121231096 GPa
- SPH25Schmitz, N. F., Ploumhans, B., & Herbst, M. F. (2025) Algorithmic differentiation for plane-wave DFT: materials design, error control and learning model parameters. arXiv:2509.07785
- Nye1985Nye, J. F. (1985). Physical Properties of Crystals. Oxford University Press. Comment: Since the elastic tensor transforms equivariantly under rotations, its numerical components depend on the chosen Cartesian coordinate frame. These tabulated patterns assume a standardized orientation of the structure with respect to conventional crystallographic axes.
- Wu2005Wu, X., Vanderbilt, D., & Hamann, D. R. (2005). Systematic treatment of displacements, strains, and electric fields in density-functional perturbation theory. Physical Review B, 72(3), 035105.