API reference
This page provides a plain list of all documented functions, structs, modules and macros in DFTK. Note that this list is neither structured, complete nor particularly clean, so it only provides rough orientation at the moment. The best reference is the code itself.
DFTK.DFTK
— ModuleDFTK –- The density-functional toolkit. Provides functionality for experimenting with plane-wave density-functional theory algorithms.
DFTK.timer
— ConstantTimerOutput object used to store DFTK timings.
DFTK.AbstractArchitecture
— TypeAbstract supertype for architectures supported by DFTK.
DFTK.AdaptiveBands
— TypeDynamically adapt number of bands to be converged to ensure that the orbitals of lowest occupation are occupied to at most occupation_threshold
. To obtain rapid convergence of the eigensolver a gap between the eigenvalues of the last occupied orbital and the last computed (but not converged) orbital of gap_min
is ensured.
DFTK.Applyχ0Model
— TypeFull χ0 application, optionally dropping terms or disabling Sternheimer. All keyword arguments passed to apply_χ0
.
DFTK.AtomicLocal
— TypeAtomic local potential defined by model.atoms
.
DFTK.AtomicNonlocal
— TypeNonlocal term coming from norm-conserving pseudopotentials in Kleinmann-Bylander form. $\text{Energy} = \sum_a \sum_{ij} \sum_{n} f_n <ψ_n|p_{ai}> D_{ij} <p_{aj}|ψ_n>.$
DFTK.BlowupAbinit
— TypeBlow-up function as used in Abinit.
DFTK.BlowupCHV
— TypeBlow-up function as proposed in https://arxiv.org/abs/2210.00442 The blow-up order of the function is fixed to ensure C^2 regularity of the energies bands away from crossings and Lipschitz continuity at crossings.
DFTK.BlowupIdentity
— TypeDefault blow-up corresponding to the standard kinetic energies.
DFTK.DielectricMixing
— TypeWe use a simplification of the Resta model DOI 10.1103/physrevb.16.2717 and set $χ_0(q) = \frac{C_0 G^2}{4π (1 - C_0 G^2 / k_{TF}^2)}$ where $C_0 = 1 - ε_r$ with $ε_r$ being the macroscopic relative permittivity. We neglect $K_\text{xc}$, such that $J^{-1} ≈ \frac{k_{TF}^2 - C_0 G^2}{ε_r k_{TF}^2 - C_0 G^2}$
By default it assumes a relative permittivity of 10 (similar to Silicon). εr == 1
is equal to SimpleMixing
and εr == Inf
to KerkerMixing
. The mixing is applied to $ρ$ and $ρ_\text{spin}$ in the same way.
DFTK.DielectricModel
— TypeA localised dielectric model for $χ_0$:
\[\sqrt{L(x)} \text{IFFT} \frac{C_0 G^2}{4π (1 - C_0 G^2 / k_{TF}^2)} \text{FFT} \sqrt{L(x)}\]
where $C_0 = 1 - ε_r$, L(r)
is a real-space localization function and otherwise the same conventions are used as in DielectricMixing
.
DFTK.DivAgradOperator
— TypeNonlocal "divAgrad" operator $-½ ∇ ⋅ (A ∇)$ where $A$ is a scalar field on the real-space grid. The $-½$ is included, such that this operator is a generalisation of the kinetic energy operator (which is obtained for $A=1$).
DFTK.ElementCohenBergstresser
— MethodElement where the interaction with electrons is modelled as in CohenBergstresser1966. Only the homonuclear lattices of the diamond structure are implemented (i.e. Si, Ge, Sn).
key
may be an element symbol (like :Si
), an atomic number (e.g. 14
) or an element name (e.g. "silicon"
)
DFTK.ElementCoulomb
— MethodElement interacting with electrons via a bare Coulomb potential (for all-electron calculations) key
may be an element symbol (like :Si
), an atomic number (e.g. 14
) or an element name (e.g. "silicon"
)
DFTK.ElementGaussian
— MethodElement interacting with electrons via a Gaussian potential. Symbol is non-mandatory.
DFTK.ElementPsp
— MethodElement interacting with electrons via a pseudopotential model. key
may be an element symbol (like :Si
), an atomic number (e.g. 14
) or an element name (e.g. "silicon"
)
DFTK.Energies
— TypeA simple struct to contain a vector of energies, and utilities to print them in a nice format.
DFTK.Entropy
— TypeEntropy term -TS, where S is the electronic entropy. Turns the energy E into the free energy F=E-TS. This is in particular useful because the free energy, not the energy, is minimized at self-consistency.
DFTK.Ewald
— TypeEwald term: electrostatic energy per unit cell of the array of point charges defined by model.atoms
in a uniform background of compensating charge yielding net neutrality.
DFTK.ExternalFromFourier
— TypeExternal potential from the (unnormalized) Fourier coefficients V(G)
G is passed in cartesian coordinates
DFTK.ExternalFromReal
— TypeExternal potential from an analytic function V
(in cartesian coordinates). No low-pass filtering is performed.
DFTK.ExternalFromValues
— TypeExternal potential given as values.
DFTK.FermiTwoStage
— TypeTwo-stage Fermi level finding algorithm starting from a Gaussian-smearing guess.
DFTK.FixedBands
— TypeIn each SCF step converge exactly n_bands_converge
, computing along the way exactly n_bands_compute
(usually a few more to ease convergence in systems with small gaps).
DFTK.FourierMultiplication
— TypeFourier space multiplication, like a kinetic energy term: (Hψ)(G) = multiplier(G) ψ(G).
DFTK.GPU
— MethodConstruct a particular GPU architecture by passing the ArrayType
DFTK.Hartree
— TypeHartree term: for a decaying potential V the energy would be
1/2 ∫ρ(x)ρ(y)V(x-y) dxdy
with the integral on x in the unit cell and of y in the whole space. For the Coulomb potential with periodic boundary conditions, this is rather
1/2 ∫ρ(x)ρ(y) G(x-y) dx dy
where G is the Green's function of the periodic Laplacian with zero mean (-Δ G = sum{R} 4π δR, integral of G zero on a unit cell).
DFTK.KerkerDosMixing
— TypeThe same as KerkerMixing
, but the Thomas-Fermi wavevector is computed from the current density of states at the Fermi level.
DFTK.KerkerMixing
— TypeKerker mixing: $J^{-1} ≈ \frac{|G|^2}{k_{TF}^2 + |G|^2}$ where $k_{TF}$ is the Thomas-Fermi wave vector. For spin-polarized calculations by default the spin density is not preconditioned. Unless a non-default value for $ΔDOS_Ω$ is specified. This value should roughly be the expected difference in density of states (per unit volume) between spin-up and spin-down.
Notes:
- Abinit calls $1/k_{TF}$ the dielectric screening length (parameter dielng)
DFTK.Kinetic
— TypeKinetic energy: 1/2 sumn fn ∫ |∇ψn|^2 * blowup(-i∇Ψ).
DFTK.Kpoint
— TypeDiscretization information for $k$-point-dependent quantities such as orbitals. More generally, a $k$-point is a block of the Hamiltonian; eg collinear spin is treated by doubling the number of kpoints.
DFTK.LazyHcat
— TypeSimple wrapper to represent a matrix formed by the concatenation of column blocks: it is mostly equivalent to hcat, but doesn't allocate the full matrix. LazyHcat only supports a few multiplication routines: furthermore, a multiplication involving this structure will always yield a plain array (and not a LazyHcat structure). LazyHcat is a lightweight subset of BlockArrays.jl's functionalities, but has the advantage to be able to store GPU Arrays (BlockArrays is heavily built on Julia's CPU Array).
DFTK.LdosModel
— TypeRepresents the LDOS-based $χ_0$ model
\[χ_0(r, r') = (-D_\text{loc}(r) δ(r, r') + D_\text{loc}(r) D_\text{loc}(r') / D)\]
where $D_\text{loc}$ is the local density of states and $D$ the density of states. For details see Herbst, Levitt 2020 arXiv:2009.01665
DFTK.LibxcDensities
— MethodCompute density in real space and its derivatives starting from ρ
DFTK.LocalNonlinearity
— TypeLocal nonlinearity, with energy ∫f(ρ) where ρ is the density
DFTK.Magnetic
— TypeMagnetic term $A⋅(-i∇)$. It is assumed (but not checked) that $∇⋅A = 0$.
DFTK.MagneticFieldOperator
— TypeMagnetic field operator A⋅(-i∇).
DFTK.Model
— MethodModel(system::AbstractSystem; kwargs...)
AtomsBase-compatible Model constructor. Sets structural information (atoms
, positions
, lattice
, n_electrons
etc.) from the passed system
.
DFTK.Model
— MethodModel(lattice, atoms, positions; n_electrons, magnetic_moments, terms, temperature,
- smearing, spin_polarization, symmetries)
Creates the physical specification of a model (without any discretization information).
n_electrons
is taken from atoms
if not specified.
spin_polarization
is :none by default (paired electrons) unless any of the elements has a non-zero initial magnetic moment. In this case the spin_polarization will be :collinear.
magnetic_moments
is only used to determine the symmetry and the spin_polarization
; it is not stored inside the datastructure.
smearing
is Fermi-Dirac if temperature
is non-zero, none otherwise
The symmetries
kwarg allows (a) to pass true
/ false
to enable / disable the automatic determination of lattice symmetries or (b) to pass an explicit list of symmetry operations to use for lowering the computational effort. The default behaviour is equal to true
, namely that the code checks the specified model in form of the Hamiltonian terms
, lattice
, atoms
and magnetic_moments
parameters and from these automatically determines a set of symmetries it can safely use. If you want to pass custom symmetry operations (e.g. a reduced or extended set) use the symmetry_operations
function. Notice that this may lead to wrong results if e.g. the external potential breaks some of the passed symmetries. Use false
to turn off symmetries completely.
DFTK.Model
— MethodModel(model; [lattice, positions, atoms, kwargs...])
-Model{T}(model; [lattice, positions, atoms, kwargs...])
Construct an identical model to model
with the option to change some of the contained parameters. This constructor is useful for changing the data type in the model or for changing lattice
or positions
in geometry/lattice optimisations.
DFTK.NbandsAlgorithm
— TypeNbandsAlgorithm subtypes determine how many bands to compute and converge in each SCF step.
DFTK.NonlocalOperator
— TypeNonlocal operator in Fourier space in Kleinman-Bylander format, defined by its projectors P matrix and coupling terms D: Hψ = PDP' ψ.
DFTK.NoopOperator
— TypeNoop operation: don't do anything. Useful for energy terms that don't depend on the orbitals at all (eg nuclei-nuclei interaction).
DFTK.PairwisePotential
— MethodPairwise terms: Pairwise potential between nuclei, e.g., Van der Waals potentials, such as Lennard—Jones terms. The potential is dependent on the distance between to atomic positions and the pairwise atomic types: For a distance d
between to atoms A
and B
, the potential is V(d, params[(A, B)])
. The parameters max_radius
is of 100
by default, and gives the maximum distance (in Cartesian coordinates) between nuclei for which we consider interactions.
DFTK.PlaneWaveBasis
— TypeA plane-wave discretized Model
. Normalization conventions:
- Things that are expressed in the G basis are normalized so that if $x$ is the vector, then the actual function is $\sum_G x_G e_G$ with $e_G(x) = e^{iG x} / \sqrt(\Omega)$, where $\Omega$ is the unit cell volume. This is so that, eg $norm(ψ) = 1$ gives the correct normalization. This also holds for the density and the potentials.
- Quantities expressed on the real-space grid are in actual values.
ifft
and fft
convert between these representations.
DFTK.PlaneWaveBasis
— MethodCreates a PlaneWaveBasis
using the kinetic energy cutoff Ecut
and a Monkhorst-Pack $k$-point grid. The MP grid can either be specified directly with kgrid
providing the number of points in each dimension and kshift
the shift (0 or 1/2 in each direction). If not specified a grid is generated using kgrid_from_minimal_spacing
with a minimal spacing of 2π * 0.022
per Bohr.
DFTK.PlaneWaveBasis
— MethodCreates a new basis identical to basis
, but with a custom set of kpoints
DFTK.PreconditionerNone
— TypeNo preconditioning
DFTK.PreconditionerTPA
— Type(simplified version of) Tetter-Payne-Allan preconditioning ↑ M.P. Teter, M.C. Payne and D.C. Allan, Phys. Rev. B 40, 12255 (1989).
DFTK.PspCorrection
— TypePseudopotential correction energy. TODO discuss the need for this.
DFTK.PspHgh
— MethodPspHgh(path[, identifier, description])
Construct a Hartwigsen, Goedecker, Teter, Hutter separable dual-space Gaussian pseudopotential (1998) from file.
DFTK.PspUpf
— MethodPspUpf(path[, identifier])
Construct a Unified Pseudopotential Format pseudopotential from file.
Does not support:
- Non-linear core correction
- Fully-realtivistic / spin-orbit pseudos
- Bare Coulomb / all-electron potentials
- Semilocal potentials
- Ultrasoft potentials
- Projector-augmented wave potentials
- GIPAW reconstruction data
DFTK.RealFourierOperator
— TypeLinear operators that act on tuples (real, fourier) The main entry point is apply!(out, op, in)
which performs the operation out += op*in where out and in are named tuples (real, fourier) They also implement mul! and Matrix(op) for exploratory use.
DFTK.RealSpaceMultiplication
— TypeReal space multiplication by a potential: (Hψ)(r) = V(r) ψ(r).
DFTK.SimpleMixing
— TypeSimple mixing: $J^{-1} ≈ 1$
DFTK.TermNoop
— TypeA term with a constant zero energy.
DFTK.Xc
— TypeExchange-correlation term, defined by a list of functionals and usually evaluated through libxc.
DFTK.χ0Mixing
— TypeGeneric mixing function using a model for the susceptibility composed of the sum of the χ0terms
. For valid χ0terms
See the subtypes of χ0Model
. The dielectric model is solved in real space using a GMRES. Either the full kernel (RPA=false
) or only the Hartree kernel (RPA=true
) are employed. verbose=true
lets the GMRES run in verbose mode (useful for debugging).
AbstractFFTs.fft!
— MethodIn-place version of fft!
. NOTE: If kpt
is given, not only f_fourier
but also f_real
is overwritten.
AbstractFFTs.fft
— Methodfft(basis::PlaneWaveBasis, [kpt::Kpoint, ] f_real)
Perform an FFT to obtain the Fourier representation of f_real
. If kpt
is given, the coefficients are truncated to the k-dependent spherical basis set.
AbstractFFTs.ifft!
— MethodIn-place version of ifft
.
AbstractFFTs.ifft
— Methodifft(basis::PlaneWaveBasis, [kpt::Kpoint, ] f_fourier)
Perform an iFFT to obtain the quantity defined by f_fourier
defined on the k-dependent spherical basis set (if kpt
is given) or the k-independent cubic (if it is not) on the real-space grid.
AtomsBase.atomic_symbol
— MethodChemical symbol corresponding to an element
AtomsBase.atomic_system
— Functionatomic_system(model::DFTK.Model, magnetic_moments=[])
-atomic_system(lattice, atoms, positions, magnetic_moments=[])
Construct an AtomsBase atomic system from a DFTK model
and associated magnetic moments or from the usual lattice
, atoms
and positions
list used in DFTK plus magnetic moments.
AtomsBase.periodic_system
— Functionperiodic_system(model::DFTK.Model, magnetic_moments=[])
-periodic_system(lattice, atoms, positions, magnetic_moments=[])
Construct an AtomsBase atomic system from a DFTK model
and associated magnetic moments or from the usual lattice
, atoms
and positions
list used in DFTK plus magnetic moments.
Brillouin.KPaths.irrfbz_path
— MethodExtract the high-symmetry $k$-point path corresponding to the passed model
using Brillouin
. Uses the conventions described in the reference work by Cracknell, Davies, Miller, and Love (CDML). Of note, this has minor differences to the $k$-path reference (Y. Himuma et. al. Comput. Mater. Sci. 128, 140 (2017)) underlying the path-choices of Brillouin.jl
, specifically for oA and mC Bravais types.
If the cell is a supercell of a smaller primitive cell, the standard $k$-path of the associated primitive cell is returned. So, the high-symmetry $k$ points are those of the primitive cell Brillouin zone, not those of the supercell Brillouin zone.
The dim
argument allows to artificially truncate the dimension of the employed model, e.g. allowing to plot a 2D bandstructure of a 3D model (useful for example for plotting band structures of sheets with dim=2
).
Due to lacking support in Spglib.jl
for two-dimensional lattices it is (a) assumed that model.lattice
is a conventional lattice and (b) required to pass the space group number using the sgnum
keyword argument.
DFTK.CROP
— FunctionCROP-accelerated root-finding iteration for f
, starting from x0
and keeping a history of m
steps. Optionally warming
specifies the number of non-accelerated steps to perform for warming up the history.
DFTK.G_vectors
— MethodG_vectors(basis::PlaneWaveBasis)
-G_vectors(basis::PlaneWaveBasis, kpt::Kpoint)
The list of wave vectors $G$ in reduced (integer) coordinates of a basis
or a $k$-point kpt
.
DFTK.G_vectors
— MethodG_vectors([architecture=AbstractArchitecture], fft_size::Tuple)
The wave vectors G
in reduced (integer) coordinates for a cubic basis set of given sizes.
DFTK.G_vectors_cart
— MethodG_vectors_cart(basis::PlaneWaveBasis)
-G_vectors_cart(basis::PlaneWaveBasis, kpt::Kpoint)
The list of $G$ vectors of a given basis
or kpt
, in cartesian coordinates.
DFTK.Gplusk_vectors
— MethodGplusk_vectors(basis::PlaneWaveBasis, kpt::Kpoint)
The list of $G + k$ vectors, in reduced coordinates.
DFTK.Gplusk_vectors_cart
— MethodGplusk_vectors_cart(basis::PlaneWaveBasis, kpt::Kpoint)
The list of $G + k$ vectors, in cartesian coordinates.
DFTK.Gplusk_vectors_in_supercell
— MethodMaps all $k+G$ vectors of an given basis as $G$ vectors of the supercell basis, in reduced coordinates.
DFTK.HybridMixing
— MethodThe model for the susceptibility is
\[\begin{aligned} +
API reference
This page provides a plain list of all documented functions, structs, modules and macros in DFTK. Note that this list is neither structured, complete nor particularly clean, so it only provides rough orientation at the moment. The best reference is the code itself.
DFTK.DFTK
— ModuleDFTK –- The density-functional toolkit. Provides functionality for experimenting with plane-wave density-functional theory algorithms.
DFTK.timer
— ConstantTimerOutput object used to store DFTK timings.
DFTK.AbstractArchitecture
— TypeAbstract supertype for architectures supported by DFTK.
DFTK.AdaptiveBands
— TypeDynamically adapt number of bands to be converged to ensure that the orbitals of lowest occupation are occupied to at most occupation_threshold
. To obtain rapid convergence of the eigensolver a gap between the eigenvalues of the last occupied orbital and the last computed (but not converged) orbital of gap_min
is ensured.
DFTK.Applyχ0Model
— TypeFull χ0 application, optionally dropping terms or disabling Sternheimer. All keyword arguments passed to apply_χ0
.
DFTK.AtomicLocal
— TypeAtomic local potential defined by model.atoms
.
DFTK.AtomicNonlocal
— TypeNonlocal term coming from norm-conserving pseudopotentials in Kleinmann-Bylander form. $\text{Energy} = \sum_a \sum_{ij} \sum_{n} f_n <ψ_n|p_{ai}> D_{ij} <p_{aj}|ψ_n>.$
DFTK.BlowupAbinit
— TypeBlow-up function as used in Abinit.
DFTK.BlowupCHV
— TypeBlow-up function as proposed in https://arxiv.org/abs/2210.00442 The blow-up order of the function is fixed to ensure C^2 regularity of the energies bands away from crossings and Lipschitz continuity at crossings.
DFTK.BlowupIdentity
— TypeDefault blow-up corresponding to the standard kinetic energies.
DFTK.DielectricMixing
— TypeWe use a simplification of the Resta model DOI 10.1103/physrevb.16.2717 and set $χ_0(q) = \frac{C_0 G^2}{4π (1 - C_0 G^2 / k_{TF}^2)}$ where $C_0 = 1 - ε_r$ with $ε_r$ being the macroscopic relative permittivity. We neglect $K_\text{xc}$, such that $J^{-1} ≈ \frac{k_{TF}^2 - C_0 G^2}{ε_r k_{TF}^2 - C_0 G^2}$
By default it assumes a relative permittivity of 10 (similar to Silicon). εr == 1
is equal to SimpleMixing
and εr == Inf
to KerkerMixing
. The mixing is applied to $ρ$ and $ρ_\text{spin}$ in the same way.
DFTK.DielectricModel
— TypeA localised dielectric model for $χ_0$:
\[\sqrt{L(x)} \text{IFFT} \frac{C_0 G^2}{4π (1 - C_0 G^2 / k_{TF}^2)} \text{FFT} \sqrt{L(x)}\]
where $C_0 = 1 - ε_r$, L(r)
is a real-space localization function and otherwise the same conventions are used as in DielectricMixing
.
DFTK.DivAgradOperator
— TypeNonlocal "divAgrad" operator $-½ ∇ ⋅ (A ∇)$ where $A$ is a scalar field on the real-space grid. The $-½$ is included, such that this operator is a generalisation of the kinetic energy operator (which is obtained for $A=1$).
DFTK.ElementCohenBergstresser
— MethodElement where the interaction with electrons is modelled as in CohenBergstresser1966. Only the homonuclear lattices of the diamond structure are implemented (i.e. Si, Ge, Sn).
key
may be an element symbol (like :Si
), an atomic number (e.g. 14
) or an element name (e.g. "silicon"
)
DFTK.ElementCoulomb
— MethodElement interacting with electrons via a bare Coulomb potential (for all-electron calculations) key
may be an element symbol (like :Si
), an atomic number (e.g. 14
) or an element name (e.g. "silicon"
)
DFTK.ElementGaussian
— MethodElement interacting with electrons via a Gaussian potential. Symbol is non-mandatory.
DFTK.ElementPsp
— MethodElement interacting with electrons via a pseudopotential model. key
may be an element symbol (like :Si
), an atomic number (e.g. 14
) or an element name (e.g. "silicon"
)
DFTK.Energies
— TypeA simple struct to contain a vector of energies, and utilities to print them in a nice format.
DFTK.Entropy
— TypeEntropy term -TS, where S is the electronic entropy. Turns the energy E into the free energy F=E-TS. This is in particular useful because the free energy, not the energy, is minimized at self-consistency.
DFTK.Ewald
— TypeEwald term: electrostatic energy per unit cell of the array of point charges defined by model.atoms
in a uniform background of compensating charge yielding net neutrality.
DFTK.ExternalFromFourier
— TypeExternal potential from the (unnormalized) Fourier coefficients V(G)
G is passed in cartesian coordinates
DFTK.ExternalFromReal
— TypeExternal potential from an analytic function V
(in cartesian coordinates). No low-pass filtering is performed.
DFTK.ExternalFromValues
— TypeExternal potential given as values.
DFTK.FermiTwoStage
— TypeTwo-stage Fermi level finding algorithm starting from a Gaussian-smearing guess.
DFTK.FixedBands
— TypeIn each SCF step converge exactly n_bands_converge
, computing along the way exactly n_bands_compute
(usually a few more to ease convergence in systems with small gaps).
DFTK.FourierMultiplication
— TypeFourier space multiplication, like a kinetic energy term: (Hψ)(G) = multiplier(G) ψ(G).
DFTK.GPU
— MethodConstruct a particular GPU architecture by passing the ArrayType
DFTK.Hartree
— TypeHartree term: for a decaying potential V the energy would be
1/2 ∫ρ(x)ρ(y)V(x-y) dxdy
with the integral on x in the unit cell and of y in the whole space. For the Coulomb potential with periodic boundary conditions, this is rather
1/2 ∫ρ(x)ρ(y) G(x-y) dx dy
where G is the Green's function of the periodic Laplacian with zero mean (-Δ G = sum{R} 4π δR, integral of G zero on a unit cell).
DFTK.KerkerDosMixing
— TypeThe same as KerkerMixing
, but the Thomas-Fermi wavevector is computed from the current density of states at the Fermi level.
DFTK.KerkerMixing
— TypeKerker mixing: $J^{-1} ≈ \frac{|G|^2}{k_{TF}^2 + |G|^2}$ where $k_{TF}$ is the Thomas-Fermi wave vector. For spin-polarized calculations by default the spin density is not preconditioned. Unless a non-default value for $ΔDOS_Ω$ is specified. This value should roughly be the expected difference in density of states (per unit volume) between spin-up and spin-down.
Notes:
- Abinit calls $1/k_{TF}$ the dielectric screening length (parameter dielng)
DFTK.Kinetic
— TypeKinetic energy: 1/2 sumn fn ∫ |∇ψn|^2 * blowup(-i∇Ψ).
DFTK.Kpoint
— TypeDiscretization information for $k$-point-dependent quantities such as orbitals. More generally, a $k$-point is a block of the Hamiltonian; eg collinear spin is treated by doubling the number of kpoints.
DFTK.LazyHcat
— TypeSimple wrapper to represent a matrix formed by the concatenation of column blocks: it is mostly equivalent to hcat, but doesn't allocate the full matrix. LazyHcat only supports a few multiplication routines: furthermore, a multiplication involving this structure will always yield a plain array (and not a LazyHcat structure). LazyHcat is a lightweight subset of BlockArrays.jl's functionalities, but has the advantage to be able to store GPU Arrays (BlockArrays is heavily built on Julia's CPU Array).
DFTK.LdosModel
— TypeRepresents the LDOS-based $χ_0$ model
\[χ_0(r, r') = (-D_\text{loc}(r) δ(r, r') + D_\text{loc}(r) D_\text{loc}(r') / D)\]
where $D_\text{loc}$ is the local density of states and $D$ the density of states. For details see Herbst, Levitt 2020 arXiv:2009.01665
DFTK.LibxcDensities
— MethodCompute density in real space and its derivatives starting from ρ
DFTK.LocalNonlinearity
— TypeLocal nonlinearity, with energy ∫f(ρ) where ρ is the density
DFTK.Magnetic
— TypeMagnetic term $A⋅(-i∇)$. It is assumed (but not checked) that $∇⋅A = 0$.
DFTK.MagneticFieldOperator
— TypeMagnetic field operator A⋅(-i∇).
DFTK.Model
— MethodModel(system::AbstractSystem; kwargs...)
AtomsBase-compatible Model constructor. Sets structural information (atoms
, positions
, lattice
, n_electrons
etc.) from the passed system
.
DFTK.Model
— MethodModel(lattice, atoms, positions; n_electrons, magnetic_moments, terms, temperature,
+ smearing, spin_polarization, symmetries)
Creates the physical specification of a model (without any discretization information).
n_electrons
is taken from atoms
if not specified.
spin_polarization
is :none by default (paired electrons) unless any of the elements has a non-zero initial magnetic moment. In this case the spin_polarization will be :collinear.
magnetic_moments
is only used to determine the symmetry and the spin_polarization
; it is not stored inside the datastructure.
smearing
is Fermi-Dirac if temperature
is non-zero, none otherwise
The symmetries
kwarg allows (a) to pass true
/ false
to enable / disable the automatic determination of lattice symmetries or (b) to pass an explicit list of symmetry operations to use for lowering the computational effort. The default behaviour is equal to true
, namely that the code checks the specified model in form of the Hamiltonian terms
, lattice
, atoms
and magnetic_moments
parameters and from these automatically determines a set of symmetries it can safely use. If you want to pass custom symmetry operations (e.g. a reduced or extended set) use the symmetry_operations
function. Notice that this may lead to wrong results if e.g. the external potential breaks some of the passed symmetries. Use false
to turn off symmetries completely.
DFTK.Model
— MethodModel(model; [lattice, positions, atoms, kwargs...])
+Model{T}(model; [lattice, positions, atoms, kwargs...])
Construct an identical model to model
with the option to change some of the contained parameters. This constructor is useful for changing the data type in the model or for changing lattice
or positions
in geometry/lattice optimisations.
DFTK.NbandsAlgorithm
— TypeNbandsAlgorithm subtypes determine how many bands to compute and converge in each SCF step.
DFTK.NonlocalOperator
— TypeNonlocal operator in Fourier space in Kleinman-Bylander format, defined by its projectors P matrix and coupling terms D: Hψ = PDP' ψ.
DFTK.NoopOperator
— TypeNoop operation: don't do anything. Useful for energy terms that don't depend on the orbitals at all (eg nuclei-nuclei interaction).
DFTK.PairwisePotential
— MethodPairwise terms: Pairwise potential between nuclei, e.g., Van der Waals potentials, such as Lennard—Jones terms. The potential is dependent on the distance between to atomic positions and the pairwise atomic types: For a distance d
between to atoms A
and B
, the potential is V(d, params[(A, B)])
. The parameters max_radius
is of 100
by default, and gives the maximum distance (in Cartesian coordinates) between nuclei for which we consider interactions.
DFTK.PlaneWaveBasis
— TypeA plane-wave discretized Model
. Normalization conventions:
- Things that are expressed in the G basis are normalized so that if $x$ is the vector, then the actual function is $\sum_G x_G e_G$ with $e_G(x) = e^{iG x} / \sqrt(\Omega)$, where $\Omega$ is the unit cell volume. This is so that, eg $norm(ψ) = 1$ gives the correct normalization. This also holds for the density and the potentials.
- Quantities expressed on the real-space grid are in actual values.
ifft
and fft
convert between these representations.
DFTK.PlaneWaveBasis
— MethodCreates a PlaneWaveBasis
using the kinetic energy cutoff Ecut
and a Monkhorst-Pack $k$-point grid. The MP grid can either be specified directly with kgrid
providing the number of points in each dimension and kshift
the shift (0 or 1/2 in each direction). If not specified a grid is generated using kgrid_from_minimal_spacing
with a minimal spacing of 2π * 0.022
per Bohr.
DFTK.PlaneWaveBasis
— MethodCreates a new basis identical to basis
, but with a custom set of kpoints
DFTK.PreconditionerNone
— TypeNo preconditioning
DFTK.PreconditionerTPA
— Type(simplified version of) Tetter-Payne-Allan preconditioning ↑ M.P. Teter, M.C. Payne and D.C. Allan, Phys. Rev. B 40, 12255 (1989).
DFTK.PspCorrection
— TypePseudopotential correction energy. TODO discuss the need for this.
DFTK.PspHgh
— MethodPspHgh(path[, identifier, description])
Construct a Hartwigsen, Goedecker, Teter, Hutter separable dual-space Gaussian pseudopotential (1998) from file.
DFTK.PspUpf
— MethodPspUpf(path[, identifier])
Construct a Unified Pseudopotential Format pseudopotential from file.
Does not support:
- Non-linear core correction
- Fully-realtivistic / spin-orbit pseudos
- Bare Coulomb / all-electron potentials
- Semilocal potentials
- Ultrasoft potentials
- Projector-augmented wave potentials
- GIPAW reconstruction data
DFTK.RealFourierOperator
— TypeLinear operators that act on tuples (real, fourier) The main entry point is apply!(out, op, in)
which performs the operation out += op*in where out and in are named tuples (real, fourier) They also implement mul! and Matrix(op) for exploratory use.
DFTK.RealSpaceMultiplication
— TypeReal space multiplication by a potential: (Hψ)(r) = V(r) ψ(r).
DFTK.SimpleMixing
— TypeSimple mixing: $J^{-1} ≈ 1$
DFTK.TermNoop
— TypeA term with a constant zero energy.
DFTK.Xc
— TypeExchange-correlation term, defined by a list of functionals and usually evaluated through libxc.
DFTK.χ0Mixing
— TypeGeneric mixing function using a model for the susceptibility composed of the sum of the χ0terms
. For valid χ0terms
See the subtypes of χ0Model
. The dielectric model is solved in real space using a GMRES. Either the full kernel (RPA=false
) or only the Hartree kernel (RPA=true
) are employed. verbose=true
lets the GMRES run in verbose mode (useful for debugging).
AbstractFFTs.fft!
— MethodIn-place version of fft!
. NOTE: If kpt
is given, not only f_fourier
but also f_real
is overwritten.
AbstractFFTs.fft
— Methodfft(basis::PlaneWaveBasis, [kpt::Kpoint, ] f_real)
Perform an FFT to obtain the Fourier representation of f_real
. If kpt
is given, the coefficients are truncated to the k-dependent spherical basis set.
AbstractFFTs.ifft!
— MethodIn-place version of ifft
.
AbstractFFTs.ifft
— Methodifft(basis::PlaneWaveBasis, [kpt::Kpoint, ] f_fourier)
Perform an iFFT to obtain the quantity defined by f_fourier
defined on the k-dependent spherical basis set (if kpt
is given) or the k-independent cubic (if it is not) on the real-space grid.
AtomsBase.atomic_symbol
— MethodChemical symbol corresponding to an element
AtomsBase.atomic_system
— Functionatomic_system(model::DFTK.Model, magnetic_moments=[])
+atomic_system(lattice, atoms, positions, magnetic_moments=[])
Construct an AtomsBase atomic system from a DFTK model
and associated magnetic moments or from the usual lattice
, atoms
and positions
list used in DFTK plus magnetic moments.
AtomsBase.periodic_system
— Functionperiodic_system(model::DFTK.Model, magnetic_moments=[])
+periodic_system(lattice, atoms, positions, magnetic_moments=[])
Construct an AtomsBase atomic system from a DFTK model
and associated magnetic moments or from the usual lattice
, atoms
and positions
list used in DFTK plus magnetic moments.
Brillouin.KPaths.irrfbz_path
— MethodExtract the high-symmetry $k$-point path corresponding to the passed model
using Brillouin
. Uses the conventions described in the reference work by Cracknell, Davies, Miller, and Love (CDML). Of note, this has minor differences to the $k$-path reference (Y. Himuma et. al. Comput. Mater. Sci. 128, 140 (2017)) underlying the path-choices of Brillouin.jl
, specifically for oA and mC Bravais types.
If the cell is a supercell of a smaller primitive cell, the standard $k$-path of the associated primitive cell is returned. So, the high-symmetry $k$ points are those of the primitive cell Brillouin zone, not those of the supercell Brillouin zone.
The dim
argument allows to artificially truncate the dimension of the employed model, e.g. allowing to plot a 2D bandstructure of a 3D model (useful for example for plotting band structures of sheets with dim=2
).
Due to lacking support in Spglib.jl
for two-dimensional lattices it is (a) assumed that model.lattice
is a conventional lattice and (b) required to pass the space group number using the sgnum
keyword argument.
DFTK.CROP
— FunctionCROP-accelerated root-finding iteration for f
, starting from x0
and keeping a history of m
steps. Optionally warming
specifies the number of non-accelerated steps to perform for warming up the history.
DFTK.G_vectors
— MethodG_vectors(basis::PlaneWaveBasis)
+G_vectors(basis::PlaneWaveBasis, kpt::Kpoint)
The list of wave vectors $G$ in reduced (integer) coordinates of a basis
or a $k$-point kpt
.
DFTK.G_vectors
— MethodG_vectors([architecture=AbstractArchitecture], fft_size::Tuple)
The wave vectors G
in reduced (integer) coordinates for a cubic basis set of given sizes.
DFTK.G_vectors_cart
— MethodG_vectors_cart(basis::PlaneWaveBasis)
+G_vectors_cart(basis::PlaneWaveBasis, kpt::Kpoint)
The list of $G$ vectors of a given basis
or kpt
, in cartesian coordinates.
DFTK.Gplusk_vectors
— MethodGplusk_vectors(basis::PlaneWaveBasis, kpt::Kpoint)
The list of $G + k$ vectors, in reduced coordinates.
DFTK.Gplusk_vectors_cart
— MethodGplusk_vectors_cart(basis::PlaneWaveBasis, kpt::Kpoint)
The list of $G + k$ vectors, in cartesian coordinates.
DFTK.Gplusk_vectors_in_supercell
— MethodMaps all $k+G$ vectors of an given basis as $G$ vectors of the supercell basis, in reduced coordinates.
DFTK.HybridMixing
— MethodThe model for the susceptibility is
\[\begin{aligned} χ_0(r, r') &= (-D_\text{loc}(r) δ(r, r') + D_\text{loc}(r) D_\text{loc}(r') / D) \\ &+ \sqrt{L(x)} \text{IFFT} \frac{C_0 G^2}{4π (1 - C_0 G^2 / k_{TF}^2)} \text{FFT} \sqrt{L(x)} -\end{aligned}\]
where $C_0 = 1 - ε_r$, $D_\text{loc}$ is the local density of states, $D$ is the density of states and the same convention for parameters are used as in DielectricMixing
. Additionally there is the real-space localization function L(r)
. For details see Herbst, Levitt 2020 arXiv:2009.01665
Important kwargs
passed on to χ0Mixing
RPA
: Is the random-phase approximation used for the kernel (i.e. only Hartree kernel is used and not XC kernel)verbose
: Run the GMRES in verbose mode.reltol
: Relative tolerance for GMRES
DFTK.IncreaseMixingTemperature
— MethodIncrease the temperature used for computing the SCF preconditioners. Initially the temperature is increased by a factor
, which is then smoothly lowered towards the temperature used within the model as the SCF converges. Once the density change is below above_ρdiff
the mixing temperature is equal to the model temperature.
DFTK.LdosMixing
— MethodThe model for the susceptibility is
\[\begin{aligned} +\end{aligned}\]
where $C_0 = 1 - ε_r$, $D_\text{loc}$ is the local density of states, $D$ is the density of states and the same convention for parameters are used as in DielectricMixing
. Additionally there is the real-space localization function L(r)
. For details see Herbst, Levitt 2020 arXiv:2009.01665
Important kwargs
passed on to χ0Mixing
RPA
: Is the random-phase approximation used for the kernel (i.e. only Hartree kernel is used and not XC kernel)verbose
: Run the GMRES in verbose mode.reltol
: Relative tolerance for GMRES
DFTK.IncreaseMixingTemperature
— MethodIncrease the temperature used for computing the SCF preconditioners. Initially the temperature is increased by a factor
, which is then smoothly lowered towards the temperature used within the model as the SCF converges. Once the density change is below above_ρdiff
the mixing temperature is equal to the model temperature.
DFTK.LdosMixing
— MethodThe model for the susceptibility is
\[\begin{aligned} χ_0(r, r') &= (-D_\text{loc}(r) δ(r, r') + D_\text{loc}(r) D_\text{loc}(r') / D) -\end{aligned}\]
where $D_\text{loc}$ is the local density of states, $D$ is the density of states. For details see Herbst, Levitt 2020 arXiv:2009.01665.
Important kwargs
passed on to χ0Mixing
RPA
: Is the random-phase approximation used for the kernel (i.e. only Hartree kernel is used and not XC kernel)verbose
: Run the GMRES in verbose mode.reltol
: Relative tolerance for GMRES
DFTK.ScfAcceptImprovingStep
— MethodAccept a step if the energy is at most increasing by max_energy
and the residual is at most max_relative_residual
times the residual in the previous step.
DFTK.ScfConvergenceDensity
— MethodFlag convergence by using the L2Norm of the change between input density and unpreconditioned output density (ρout)
DFTK.ScfConvergenceEnergy
— MethodFlag convergence as soon as total energy change drops below tolerance
DFTK.ScfConvergenceForce
— MethodFlag convergence on the change in cartesian force between two iterations.
DFTK.ScfDefaultCallback
— MethodDefault callback function for self_consistent_field
and newton
, which prints a convergence table.
DFTK.ScfDiagtol
— MethodDetermine the tolerance used for the next diagonalization. This function takes $|ρnext - ρin|$ and multiplies it with ratio_ρdiff
to get the next diagtol
, ensuring additionally that the returned value is between diagtol_min
and diagtol_max
and never increases.
DFTK.ScfPlotTrace
— FunctionPlot the trace of an SCF, i.e. the absolute error of the total energy at each iteration versus the converged energy in a semilog plot. By default a new plot canvas is generated, but an existing one can be passed and reused along with kwargs
for the call to plot!
. Requires Plots to be loaded.
DFTK.ScfSaveCheckpoints
— FunctionAdds simplistic checkpointing to a DFTK self-consistent field calculation. Requires JLD2 to be loaded.
DFTK.apply_K
— Methodapply_K(basis::PlaneWaveBasis, δψ, ψ, ρ, occupation)
Compute the application of K defined at ψ to δψ. ρ is the density issued from ψ. δψ also generates a δρ, computed with compute_δρ
.
DFTK.apply_kernel
— Methodapply_kernel(basis::PlaneWaveBasis, δρ; kwargs...)
Computes the potential response to a perturbation δρ in real space, as a 4D (i,j,k,σ) array.
DFTK.apply_symop
— MethodApply a symmetry operation to eigenvectors ψk
at a given kpoint
to obtain an equivalent point in [-0.5, 0.5)^3 and associated eigenvectors (expressed in the basis of the new $k$-point).
DFTK.apply_symop
— MethodApply a symmetry operation to a density.
DFTK.apply_Ω
— Methodapply_Ω(δψ, ψ, H::Hamiltonian, Λ)
Compute the application of Ω defined at ψ to δψ. H is the Hamiltonian computed from ψ and Λ is the set of Rayleigh coefficients ψk' * Hk * ψk at each k-point.
DFTK.apply_χ0
— MethodGet the density variation δρ corresponding to a potential variation δV.
DFTK.attach_psp
— Methodattach_psp(system::AbstractSystem, pspmap::AbstractDict{Symbol,String})
-attach_psp(system::AbstractSystem; psps::String...)
Return a new system with the pseudopotential
property of all atoms set according to the passed pspmap
, which maps from the atomic symbol to a pseudopotential identifier. Alternatively the mapping from atomic symbol to pseudopotential identifier can also be passed as keyword arguments. An empty string can be used to denote elements where the full Coulomb potential should be employed.
Examples
Select pseudopotentials for all silicon and oxygen atoms in the system.
julia> attach_psp(system, Dict(:Si => "hgh/lda/si-q4", :O => "hgh/lda/o-q6")
Same thing but using the kwargs syntax:
julia> attach_psp(system, Si="hgh/lda/si-q4", O="hgh/lda/o-q6")
DFTK.build_fft_plans!
— MethodPlan a FFT of type T
and size fft_size
, spending some time on finding an optimal algorithm. (Inplace, out-of-place) x (forward, backward) FFT plans are returned.
DFTK.build_form_factors
— MethodBuild form factors (Fourier transforms of projectors) for an atom centered at 0.
DFTK.build_projection_vectors_
— MethodBuild projection vectors for a atoms array generated by term_nonlocal
\[\begin{aligned} +\end{aligned}\]
where $D_\text{loc}$ is the local density of states, $D$ is the density of states. For details see Herbst, Levitt 2020 arXiv:2009.01665.
Important kwargs
passed on to χ0Mixing
RPA
: Is the random-phase approximation used for the kernel (i.e. only Hartree kernel is used and not XC kernel)verbose
: Run the GMRES in verbose mode.reltol
: Relative tolerance for GMRES
DFTK.ScfAcceptImprovingStep
— MethodAccept a step if the energy is at most increasing by max_energy
and the residual is at most max_relative_residual
times the residual in the previous step.
DFTK.ScfConvergenceDensity
— MethodFlag convergence by using the L2Norm of the change between input density and unpreconditioned output density (ρout)
DFTK.ScfConvergenceEnergy
— MethodFlag convergence as soon as total energy change drops below tolerance
DFTK.ScfConvergenceForce
— MethodFlag convergence on the change in cartesian force between two iterations.
DFTK.ScfDefaultCallback
— MethodDefault callback function for self_consistent_field
and newton
, which prints a convergence table.
DFTK.ScfDiagtol
— MethodDetermine the tolerance used for the next diagonalization. This function takes $|ρnext - ρin|$ and multiplies it with ratio_ρdiff
to get the next diagtol
, ensuring additionally that the returned value is between diagtol_min
and diagtol_max
and never increases.
DFTK.ScfPlotTrace
— FunctionPlot the trace of an SCF, i.e. the absolute error of the total energy at each iteration versus the converged energy in a semilog plot. By default a new plot canvas is generated, but an existing one can be passed and reused along with kwargs
for the call to plot!
. Requires Plots to be loaded.
DFTK.ScfSaveCheckpoints
— FunctionAdds simplistic checkpointing to a DFTK self-consistent field calculation. Requires JLD2 to be loaded.
DFTK.apply_K
— Methodapply_K(basis::PlaneWaveBasis, δψ, ψ, ρ, occupation)
Compute the application of K defined at ψ to δψ. ρ is the density issued from ψ. δψ also generates a δρ, computed with compute_δρ
.
DFTK.apply_kernel
— Methodapply_kernel(basis::PlaneWaveBasis, δρ; kwargs...)
Computes the potential response to a perturbation δρ in real space, as a 4D (i,j,k,σ) array.
DFTK.apply_symop
— MethodApply a symmetry operation to eigenvectors ψk
at a given kpoint
to obtain an equivalent point in [-0.5, 0.5)^3 and associated eigenvectors (expressed in the basis of the new $k$-point).
DFTK.apply_symop
— MethodApply a symmetry operation to a density.
DFTK.apply_Ω
— Methodapply_Ω(δψ, ψ, H::Hamiltonian, Λ)
Compute the application of Ω defined at ψ to δψ. H is the Hamiltonian computed from ψ and Λ is the set of Rayleigh coefficients ψk' * Hk * ψk at each k-point.
DFTK.apply_χ0
— MethodGet the density variation δρ corresponding to a potential variation δV.
DFTK.attach_psp
— Methodattach_psp(system::AbstractSystem, pspmap::AbstractDict{Symbol,String})
+attach_psp(system::AbstractSystem; psps::String...)
Return a new system with the pseudopotential
property of all atoms set according to the passed pspmap
, which maps from the atomic symbol to a pseudopotential identifier. Alternatively the mapping from atomic symbol to pseudopotential identifier can also be passed as keyword arguments. An empty string can be used to denote elements where the full Coulomb potential should be employed.
Examples
Select pseudopotentials for all silicon and oxygen atoms in the system.
julia> attach_psp(system, Dict(:Si => "hgh/lda/si-q4", :O => "hgh/lda/o-q6")
Same thing but using the kwargs syntax:
julia> attach_psp(system, Si="hgh/lda/si-q4", O="hgh/lda/o-q6")
DFTK.build_fft_plans!
— MethodPlan a FFT of type T
and size fft_size
, spending some time on finding an optimal algorithm. (Inplace, out-of-place) x (forward, backward) FFT plans are returned.
DFTK.build_form_factors
— MethodBuild form factors (Fourier transforms of projectors) for an atom centered at 0.
DFTK.build_projection_vectors_
— MethodBuild projection vectors for a atoms array generated by term_nonlocal
\[\begin{aligned} H_{\rm at} &= \sum_{ij} C_{ij} \ket{p_i} \bra{p_j} \\ H_{\rm per} &= \sum_R \sum_{ij} C_{ij} \ket{p_i(x-R)} \bra{p_j(x-R)} \end{aligned}\]
\[\begin{aligned} \braket{e_k(G') \middle| H_{\rm per}}{e_k(G)} &= \ldots \\ &= \frac{1}{Ω} \sum_{ij} C_{ij} \hat p_i(k+G') \hat p_j^*(k+G), -\end{aligned}\]
where $\hat p_i(q) = ∫_{ℝ^3} p_i(r) e^{-iq·r} dr$.
We store $\frac{1}{\sqrt Ω} \hat p_i(k+G)$ in proj_vectors
.
DFTK.bzmesh_ir_wedge
— Method bzmesh_ir_wedge(kgrid_size, symmetries; kshift=[0, 0, 0])
Construct the irreducible wedge of a uniform Brillouin zone mesh for sampling $k$-points, given the crystal symmetries symmetries
. Returns the list of irreducible $k$-point (fractional) coordinates, the associated weights and the new symmetries
compatible with the grid.
DFTK.bzmesh_uniform
— Methodbzmesh_uniform(kgrid_size; kshift=[0, 0, 0])
Construct a (shifted) uniform Brillouin zone mesh for sampling the $k$-points. Returns all $k$-point coordinates, appropriate weights and the identity SymOp.
DFTK.cell_to_supercell
— MethodTranspose all data from a given self-consistent-field result from unit cell to supercell conventions. The parameters to adapt are the following:
basis_supercell
andψ_supercell
are computed by the routines above.- The supercell occupations vector is the concatenation of all input occupations vectors.
- The supercell density is computed with supercell occupations and
ψ_supercell
. - Supercell energies are the multiplication of input energies by the number of unit cells in the supercell.
Other parameters stay untouched.
DFTK.cell_to_supercell
— MethodConstruct a plane-wave basis whose unit cell is the supercell associated to an input basis $k$-grid. All other parameters are modified so that the respective physical systems associated to both basis are equivalent.
DFTK.cell_to_supercell
— MethodRe-organize Bloch waves computed in a given basis as Bloch waves of the associated supercell basis. The output ψ_supercell
have a single component at $Γ$-point, such that ψ_supercell[Γ][:, k+n]
contains ψ[k][:, n]
, within normalization on the supercell.
DFTK.cg!
— MethodImplementation of the conjugate gradient method which allows for preconditioning and projection operations along iterations.
DFTK.charge_ionic
— MethodReturn the total ionic charge of an atom type (nuclear charge - core electrons)
DFTK.charge_nuclear
— MethodReturn the total nuclear charge of an atom type
DFTK.cis2pi
— MethodFunction to compute exp(2π i x)
DFTK.compute_Ak_gaussian_guess
— MethodCompute the matrix $[A_k]_{m,n} = \langle ψ_m^k | g^{\text{per}}_n \rangle$
$g^{per}_n$ are periodized gaussians whose respective centers are given as an (num_bands,1) array [ [center 1], ... ].
Centers are to be given in lattice coordinates and G_vectors in reduced coordinates. The dot product is computed in the Fourier space.
Given an orbital $g_n$, the periodized orbital is defined by : $g^{per}_n = \sum\limits_{R \in {\rm lattice}} g_n( \cdot - R)$. The Fourier coefficient of $g^{per}_n$ at any G is given by the value of the Fourier transform of $g_n$ in G.
DFTK.compute_current
— MethodComputes the probability (not charge) current, ∑ fn Im(ψn* ∇ψn)
DFTK.compute_density
— Methodcompute_density(basis::PlaneWaveBasis, ψ::AbstractVector, occupation::AbstractVector)
Compute the density for a wave function ψ
discretized on the plane-wave grid basis
, where the individual k-points are occupied according to occupation
. ψ
should be one coefficient matrix per $k$-point. It is possible to ask only for occupations higher than a certain level to be computed by using an optional occupation_threshold
. By default all occupation numbers are considered.
DFTK.compute_dos
— MethodTotal density of states at energy ε
DFTK.compute_dynmat
— MethodCompute the dynamical matrix in the form of a $3×n_{ m atoms}×3×n_{ m atoms}$ tensor in reduced coordinates.
DFTK.compute_dynmat_cart
— MethodCartesian form of compute_dynmat
.
DFTK.compute_fft_size
— MethodDetermine the minimal grid size for the cubic basis set to be able to represent product of orbitals (with the default supersampling=2
).
Optionally optimize the grid afterwards for the FFT procedure by ensuring factorization into small primes.
The function will determine the smallest parallelepiped containing the wave vectors $|G|^2/2 \leq E_\text{cut} ⋅ \text{supersampling}^2$. For an exact representation of the density resulting from wave functions represented in the spherical basis sets, supersampling
should be at least 2
.
If factors
is not empty, ensure that the resulting fft_size contains all the factors
DFTK.compute_forces
— MethodCompute the forces of an obtained SCF solution. Returns the forces wrt. the fractional lattice vectors. To get cartesian forces use compute_forces_cart
. Returns a list of lists of forces (as SVector{3}) in the same order as the atoms
and positions
in the underlying Model
.
DFTK.compute_forces_cart
— MethodCompute the cartesian forces of an obtained SCF solution in Hartree / Bohr. Returns a list of lists of forces [[force for atom in positions] for (element, positions) in atoms]
which has the same structure as the atoms
object passed to the underlying Model
.
DFTK.compute_inverse_lattice
— MethodCompute the inverse of the lattice. Takes special care of 1D or 2D cases.
DFTK.compute_kernel
— Methodcompute_kernel(basis::PlaneWaveBasis; kwargs...)
Computes a matrix representation of the full response kernel (derivative of potential with respect to density) in real space. For non-spin-polarized calculations the matrix dimension is prod(basis.fft_size)
× prod(basis.fft_size)
and for collinear spin-polarized cases it is 2prod(basis.fft_size)
× 2prod(basis.fft_size)
. In this case the matrix has effectively 4 blocks
\[\left(\begin{array}{cc} +\end{aligned}\]
where $\hat p_i(q) = ∫_{ℝ^3} p_i(r) e^{-iq·r} dr$.
We store $\frac{1}{\sqrt Ω} \hat p_i(k+G)$ in proj_vectors
.
DFTK.bzmesh_ir_wedge
— Method bzmesh_ir_wedge(kgrid_size, symmetries; kshift=[0, 0, 0])
Construct the irreducible wedge of a uniform Brillouin zone mesh for sampling $k$-points, given the crystal symmetries symmetries
. Returns the list of irreducible $k$-point (fractional) coordinates, the associated weights and the new symmetries
compatible with the grid.
DFTK.bzmesh_uniform
— Methodbzmesh_uniform(kgrid_size; kshift=[0, 0, 0])
Construct a (shifted) uniform Brillouin zone mesh for sampling the $k$-points. Returns all $k$-point coordinates, appropriate weights and the identity SymOp.
DFTK.cell_to_supercell
— MethodTranspose all data from a given self-consistent-field result from unit cell to supercell conventions. The parameters to adapt are the following:
basis_supercell
andψ_supercell
are computed by the routines above.- The supercell occupations vector is the concatenation of all input occupations vectors.
- The supercell density is computed with supercell occupations and
ψ_supercell
. - Supercell energies are the multiplication of input energies by the number of unit cells in the supercell.
Other parameters stay untouched.
DFTK.cell_to_supercell
— MethodConstruct a plane-wave basis whose unit cell is the supercell associated to an input basis $k$-grid. All other parameters are modified so that the respective physical systems associated to both basis are equivalent.
DFTK.cell_to_supercell
— MethodRe-organize Bloch waves computed in a given basis as Bloch waves of the associated supercell basis. The output ψ_supercell
have a single component at $Γ$-point, such that ψ_supercell[Γ][:, k+n]
contains ψ[k][:, n]
, within normalization on the supercell.
DFTK.cg!
— MethodImplementation of the conjugate gradient method which allows for preconditioning and projection operations along iterations.
DFTK.charge_ionic
— MethodReturn the total ionic charge of an atom type (nuclear charge - core electrons)
DFTK.charge_nuclear
— MethodReturn the total nuclear charge of an atom type
DFTK.cis2pi
— MethodFunction to compute exp(2π i x)
DFTK.compute_Ak_gaussian_guess
— MethodCompute the matrix $[A_k]_{m,n} = \langle ψ_m^k | g^{\text{per}}_n \rangle$
$g^{per}_n$ are periodized gaussians whose respective centers are given as an (num_bands,1) array [ [center 1], ... ].
Centers are to be given in lattice coordinates and G_vectors in reduced coordinates. The dot product is computed in the Fourier space.
Given an orbital $g_n$, the periodized orbital is defined by : $g^{per}_n = \sum\limits_{R \in {\rm lattice}} g_n( \cdot - R)$. The Fourier coefficient of $g^{per}_n$ at any G is given by the value of the Fourier transform of $g_n$ in G.
DFTK.compute_current
— MethodComputes the probability (not charge) current, ∑ fn Im(ψn* ∇ψn)
DFTK.compute_density
— Methodcompute_density(basis::PlaneWaveBasis, ψ::AbstractVector, occupation::AbstractVector)
Compute the density for a wave function ψ
discretized on the plane-wave grid basis
, where the individual k-points are occupied according to occupation
. ψ
should be one coefficient matrix per $k$-point. It is possible to ask only for occupations higher than a certain level to be computed by using an optional occupation_threshold
. By default all occupation numbers are considered.
DFTK.compute_dos
— MethodTotal density of states at energy ε
DFTK.compute_dynmat
— MethodCompute the dynamical matrix in the form of a $3×n_{ m atoms}×3×n_{ m atoms}$ tensor in reduced coordinates.
DFTK.compute_dynmat_cart
— MethodCartesian form of compute_dynmat
.
DFTK.compute_fft_size
— MethodDetermine the minimal grid size for the cubic basis set to be able to represent product of orbitals (with the default supersampling=2
).
Optionally optimize the grid afterwards for the FFT procedure by ensuring factorization into small primes.
The function will determine the smallest parallelepiped containing the wave vectors $|G|^2/2 \leq E_\text{cut} ⋅ \text{supersampling}^2$. For an exact representation of the density resulting from wave functions represented in the spherical basis sets, supersampling
should be at least 2
.
If factors
is not empty, ensure that the resulting fft_size contains all the factors
DFTK.compute_forces
— MethodCompute the forces of an obtained SCF solution. Returns the forces wrt. the fractional lattice vectors. To get cartesian forces use compute_forces_cart
. Returns a list of lists of forces (as SVector{3}) in the same order as the atoms
and positions
in the underlying Model
.
DFTK.compute_forces_cart
— MethodCompute the cartesian forces of an obtained SCF solution in Hartree / Bohr. Returns a list of lists of forces [[force for atom in positions] for (element, positions) in atoms]
which has the same structure as the atoms
object passed to the underlying Model
.
DFTK.compute_inverse_lattice
— MethodCompute the inverse of the lattice. Takes special care of 1D or 2D cases.
DFTK.compute_kernel
— Methodcompute_kernel(basis::PlaneWaveBasis; kwargs...)
Computes a matrix representation of the full response kernel (derivative of potential with respect to density) in real space. For non-spin-polarized calculations the matrix dimension is prod(basis.fft_size)
× prod(basis.fft_size)
and for collinear spin-polarized cases it is 2prod(basis.fft_size)
× 2prod(basis.fft_size)
. In this case the matrix has effectively 4 blocks
\[\left(\begin{array}{cc} K_{αα} & K_{αβ}\\ K_{βα} & K_{ββ} -\end{array}\right)\]
DFTK.compute_ldos
— MethodLocal density of states, in real space. weight_threshold
is a threshold to screen away small contributions to the LDOS.
DFTK.compute_occupation
— MethodCompute occupation given eigenvalues and Fermi level
DFTK.compute_occupation
— MethodCompute occupation and Fermi level given eigenvalues and using fermialg
. The tol_n_elec
gives the accuracy on the electron count which should be at least achieved.
DFTK.compute_recip_lattice
— MethodCompute the reciprocal lattice. We use the convention that the reciprocal lattice is the set of G vectors such that G ⋅ R ∈ 2π ℤ for all R in the lattice.
DFTK.compute_stresses_cart
— MethodCompute the stresses (= 1/Vol dE/d(M*lattice), taken at M=I) of an obtained SCF solution.
DFTK.compute_transfer_matrix
— MethodReturn a sparse matrix that maps quantities given on basis_in
and kpt_in
to quantities on basis_out
and kpt_out
.
DFTK.compute_transfer_matrix
— MethodReturn a list of sparse matrices (one per $k$-point) that map quantities given in the basis_in
basis to quantities given in the basis_out
basis.
DFTK.compute_unit_cell_volume
— MethodCompute unit cell volume volume. In case of 1D or 2D case, the volume is the length/surface.
DFTK.compute_δocc!
— MethodCompute the derivatives of the occupations (and of the Fermi level). The derivatives of the occupations are in-place stored in δocc. The tuple (; δocc, δεF) is returned. It is assumed the passed δocc
are initialised to zero.
DFTK.compute_δψ!
— MethodPerform in-place computations of the derivatives of the wave functions by solving a Sternheimer equation for each k
-points. It is assumed the passed δψ
are initialised to zero.
DFTK.compute_χ0
— MethodCompute the independent-particle susceptibility. Will blow up for large systems. For non-spin-polarized calculations the matrix dimension is prod(basis.fft_size)
× prod(basis.fft_size)
and for collinear spin-polarized cases it is 2prod(basis.fft_size)
× 2prod(basis.fft_size)
. In this case the matrix has effectively 4 blocks, which are:
\[\left(\begin{array}{cc} +\end{array}\right)\]
DFTK.compute_ldos
— MethodLocal density of states, in real space. weight_threshold
is a threshold to screen away small contributions to the LDOS.
DFTK.compute_occupation
— MethodCompute occupation given eigenvalues and Fermi level
DFTK.compute_occupation
— MethodCompute occupation and Fermi level given eigenvalues and using fermialg
. The tol_n_elec
gives the accuracy on the electron count which should be at least achieved.
DFTK.compute_recip_lattice
— MethodCompute the reciprocal lattice. We use the convention that the reciprocal lattice is the set of G vectors such that G ⋅ R ∈ 2π ℤ for all R in the lattice.
DFTK.compute_stresses_cart
— MethodCompute the stresses (= 1/Vol dE/d(M*lattice), taken at M=I) of an obtained SCF solution.
DFTK.compute_transfer_matrix
— MethodReturn a sparse matrix that maps quantities given on basis_in
and kpt_in
to quantities on basis_out
and kpt_out
.
DFTK.compute_transfer_matrix
— MethodReturn a list of sparse matrices (one per $k$-point) that map quantities given in the basis_in
basis to quantities given in the basis_out
basis.
DFTK.compute_unit_cell_volume
— MethodCompute unit cell volume volume. In case of 1D or 2D case, the volume is the length/surface.
DFTK.compute_δocc!
— MethodCompute the derivatives of the occupations (and of the Fermi level). The derivatives of the occupations are in-place stored in δocc. The tuple (; δocc, δεF) is returned. It is assumed the passed δocc
are initialised to zero.
DFTK.compute_δψ!
— MethodPerform in-place computations of the derivatives of the wave functions by solving a Sternheimer equation for each k
-points. It is assumed the passed δψ
are initialised to zero.
DFTK.compute_χ0
— MethodCompute the independent-particle susceptibility. Will blow up for large systems. For non-spin-polarized calculations the matrix dimension is prod(basis.fft_size)
× prod(basis.fft_size)
and for collinear spin-polarized cases it is 2prod(basis.fft_size)
× 2prod(basis.fft_size)
. In this case the matrix has effectively 4 blocks, which are:
\[\left(\begin{array}{cc} (χ_0)_{αα} & (χ_0)_{αβ} \\ (χ_0)_{βα} & (χ_0)_{ββ} -\end{array}\right)\]
DFTK.cos2pi
— MethodFunction to compute cos(2π x)
DFTK.count_n_proj
— Methodcount_n_proj(psps, psp_positions)
Number of projector functions for all angular momenta up to psp.lmax
and for all atoms in the system, including angular parts from -m:m
.
DFTK.count_n_proj
— Methodcount_n_proj(psp, l)
Number of projector functions for angular momentum l
, including angular parts from -m:m
.
DFTK.count_n_proj
— Methodcount_n_proj(psp)
Number of projector functions for all angular momenta up to psp.lmax
, including angular parts from -m:m
.
DFTK.count_n_proj_radial
— Methodcount_n_proj_radial(psp, l)
Number of projector radial functions at angular momentum l
.
DFTK.count_n_proj_radial
— Methodcount_n_proj_radial(psp)
Number of projector radial functions at all angular momenta up to psp.lmax
.
DFTK.create_supercell
— MethodConstruct a supercell of size supercell_size
from a unit cell described by its lattice
, atoms
and their positions
.
DFTK.datadir_psp
— MethodReturn the data directory with pseudopotential files
DFTK.default_fermialg
— MethodDefault selection of a Fermi level determination algorithm
DFTK.default_symmetries
— MethodDefault logic to determine the symmetry operations to be used in the model.
DFTK.default_wannier_centres
— MethodDefault random Gaussian guess for maximally-localised wannier functions generated in reduced coordinates.
DFTK.determine_spin_polarization
— Method:none
if no element has a magnetic moment, else :collinear
or :full
.
DFTK.diagonalize_all_kblocks
— MethodFunction for diagonalising each $k$-Point blow of ham one step at a time. Some logic for interpolating between $k$-points is used if interpolate_kpoints
is true and if no guesses are given. eigensolver
is the iterative eigensolver that really does the work, operating on a single $k$-Block. eigensolver
should support the API eigensolver(A, X0; prec, tol, maxiter)
prec_type
should be a function that returns a preconditioner when called as prec(ham, kpt)
DFTK.diameter
— MethodCompute the diameter of the unit cell
DFTK.direct_minimization
— MethodComputes the ground state by direct minimization. kwargs...
are passed to Optim.Options()
. Note that the resulting ψ are not necessarily eigenvectors of the Hamiltonian.
DFTK.disable_threading
— MethodConvenience function to disable all threading in DFTK and assert that Julia threading is off as well.
DFTK.divergence_real
— MethodCompute divergence of an operand function, which returns the cartesian x,y,z components in real space when called with the arguments 1 to 3. The divergence is also returned as a real-space array.
DFTK.energy_forces_ewald
— MethodStandard computation of energy and forces.
DFTK.energy_forces_ewald
— MethodComputation for phonons; required to build the dynamical matrix.
DFTK.energy_forces_ewald
— MethodCompute the electrostatic energy and forces. The energy is the electrostatic interaction energy per unit cell between point charges in a uniform background of compensating charge to yield net neutrality. The forces is the opposite of the derivative of the energy with respect to positions
.
lattice
should contain the lattice vectors as columns. charges
and positions
are the point charges and their positions (as an array of arrays) in fractional coordinates.
For now this function returns zero energy and force on non-3D systems. Use a pairwise potential term if you want to customise this treatment.
DFTK.energy_forces_pairwise
— MethodStandard computation of energy and forces.
DFTK.energy_forces_pairwise
— MethodComputation for phonons; required to build the dynamical matrix.
DFTK.energy_forces_pairwise
— MethodCompute the pairwise energy and forces. The energy is the interaction energy per unit cell between atomic sites. The forces is the opposite of the derivative of the energy with respect to positions
.
lattice
should contain the lattice vectors as columns. symbols
and positions
are the atomic elements and their positions (as an array of arrays) in fractional coordinates. V
and params
are the pairwise potential and its set of parameters (that depends on pairs of symbols).
The potential is expected to decrease quickly at infinity.
DFTK.energy_psp_correction
— MethodCompute the correction term for properly modelling the interaction of the pseudopotential core with the compensating background charge induced by the Ewald
term.
DFTK.enforce_real!
— MethodEnsure its real-space equivalent of passed Fourier-space representation is entirely real by removing wavevectors G
that don't have a -G
counterpart in the basis.
DFTK.estimate_integer_lattice_bounds
— MethodEstimate integer bounds for dense space loops from a given inequality ||Mx|| ≤ δ. For 1D and 2D systems the limit will be zero in the auxiliary dimensions.
DFTK.eval_psp_density_core_fourier
— Methodeval_psp_density_core_fourier(psp, q)
Evaluate the atomic core charge density in reciprocal space: ρval(q) = ∫{R^3} ρcore(r) e^{-iqr} dr = 4π ∫{R+} ρcore(r) sin(qr)/qr r^2 dr
DFTK.eval_psp_density_core_real
— Methodeval_psp_density_core_real(psp, r)
Evaluate the atomic core charge density in real space.
DFTK.eval_psp_density_valence_fourier
— Methodeval_psp_density_valence_fourier(psp, q)
Evaluate the atomic valence charge density in reciprocal space: ρval(q) = ∫{R^3} ρval(r) e^{-iqr} dr = 4π ∫{R+} ρval(r) sin(qr)/qr r^2 dr
DFTK.eval_psp_density_valence_real
— Methodeval_psp_density_valence_real(psp, r)
Evaluate the atomic valence charge density in real space.
DFTK.eval_psp_energy_correction
— Functioneval_psp_energy_correction([T=Float64,] psp, n_electrons)
Evaluate the energy correction to the Ewald electrostatic interaction energy of one unit cell, which is required compared the Ewald expression for point-like nuclei. n_electrons
is the number of electrons per unit cell. This defines the uniform compensating background charge, which is assumed here.
Notice: The returned result is the energy per unit cell and not the energy per volume. To obtain the latter, the caller needs to divide by the unit cell volume.
The energy correction is defined as the limit of the Fourier-transform of the local potential as $q \to 0$, using the same correction as in the Fourier-transform of the local potential: math \lim_{q \to 0} 4π N_{\rm elec} ∫_{ℝ_+} (V(r) - C(r)) \frac{\sin(qr)}{qr} r^2 dr + F[C(r)] = 4π N_{\rm elec} ∫_{ℝ_+} (V(r) + Z/r) r^2 dr
DFTK.eval_psp_local_fourier
— Methodeval_psp_local_fourier(psp, q)
Evaluate the local part of the pseudopotential in reciprocal space:
\[\begin{aligned} +\end{array}\right)\]
DFTK.cos2pi
— MethodFunction to compute cos(2π x)
DFTK.count_n_proj
— Methodcount_n_proj(psps, psp_positions)
Number of projector functions for all angular momenta up to psp.lmax
and for all atoms in the system, including angular parts from -m:m
.
DFTK.count_n_proj
— Methodcount_n_proj(psp, l)
Number of projector functions for angular momentum l
, including angular parts from -m:m
.
DFTK.count_n_proj
— Methodcount_n_proj(psp)
Number of projector functions for all angular momenta up to psp.lmax
, including angular parts from -m:m
.
DFTK.count_n_proj_radial
— Methodcount_n_proj_radial(psp, l)
Number of projector radial functions at angular momentum l
.
DFTK.count_n_proj_radial
— Methodcount_n_proj_radial(psp)
Number of projector radial functions at all angular momenta up to psp.lmax
.
DFTK.create_supercell
— MethodConstruct a supercell of size supercell_size
from a unit cell described by its lattice
, atoms
and their positions
.
DFTK.datadir_psp
— MethodReturn the data directory with pseudopotential files
DFTK.default_fermialg
— MethodDefault selection of a Fermi level determination algorithm
DFTK.default_symmetries
— MethodDefault logic to determine the symmetry operations to be used in the model.
DFTK.default_wannier_centres
— MethodDefault random Gaussian guess for maximally-localised wannier functions generated in reduced coordinates.
DFTK.determine_spin_polarization
— Method:none
if no element has a magnetic moment, else :collinear
or :full
.
DFTK.diagonalize_all_kblocks
— MethodFunction for diagonalising each $k$-Point blow of ham one step at a time. Some logic for interpolating between $k$-points is used if interpolate_kpoints
is true and if no guesses are given. eigensolver
is the iterative eigensolver that really does the work, operating on a single $k$-Block. eigensolver
should support the API eigensolver(A, X0; prec, tol, maxiter)
prec_type
should be a function that returns a preconditioner when called as prec(ham, kpt)
DFTK.diameter
— MethodCompute the diameter of the unit cell
DFTK.direct_minimization
— MethodComputes the ground state by direct minimization. kwargs...
are passed to Optim.Options()
. Note that the resulting ψ are not necessarily eigenvectors of the Hamiltonian.
DFTK.disable_threading
— MethodConvenience function to disable all threading in DFTK and assert that Julia threading is off as well.
DFTK.divergence_real
— MethodCompute divergence of an operand function, which returns the cartesian x,y,z components in real space when called with the arguments 1 to 3. The divergence is also returned as a real-space array.
DFTK.energy_forces_ewald
— MethodStandard computation of energy and forces.
DFTK.energy_forces_ewald
— MethodComputation for phonons; required to build the dynamical matrix.
DFTK.energy_forces_ewald
— MethodCompute the electrostatic energy and forces. The energy is the electrostatic interaction energy per unit cell between point charges in a uniform background of compensating charge to yield net neutrality. The forces is the opposite of the derivative of the energy with respect to positions
.
lattice
should contain the lattice vectors as columns. charges
and positions
are the point charges and their positions (as an array of arrays) in fractional coordinates.
For now this function returns zero energy and force on non-3D systems. Use a pairwise potential term if you want to customise this treatment.
DFTK.energy_forces_pairwise
— MethodStandard computation of energy and forces.
DFTK.energy_forces_pairwise
— MethodComputation for phonons; required to build the dynamical matrix.
DFTK.energy_forces_pairwise
— MethodCompute the pairwise energy and forces. The energy is the interaction energy per unit cell between atomic sites. The forces is the opposite of the derivative of the energy with respect to positions
.
lattice
should contain the lattice vectors as columns. symbols
and positions
are the atomic elements and their positions (as an array of arrays) in fractional coordinates. V
and params
are the pairwise potential and its set of parameters (that depends on pairs of symbols).
The potential is expected to decrease quickly at infinity.
DFTK.energy_psp_correction
— MethodCompute the correction term for properly modelling the interaction of the pseudopotential core with the compensating background charge induced by the Ewald
term.
DFTK.enforce_real!
— MethodEnsure its real-space equivalent of passed Fourier-space representation is entirely real by removing wavevectors G
that don't have a -G
counterpart in the basis.
DFTK.estimate_integer_lattice_bounds
— MethodEstimate integer bounds for dense space loops from a given inequality ||Mx|| ≤ δ. For 1D and 2D systems the limit will be zero in the auxiliary dimensions.
DFTK.eval_psp_density_core_fourier
— Methodeval_psp_density_core_fourier(psp, q)
Evaluate the atomic core charge density in reciprocal space: ρval(q) = ∫{R^3} ρcore(r) e^{-iqr} dr = 4π ∫{R+} ρcore(r) sin(qr)/qr r^2 dr
DFTK.eval_psp_density_core_real
— Methodeval_psp_density_core_real(psp, r)
Evaluate the atomic core charge density in real space.
DFTK.eval_psp_density_valence_fourier
— Methodeval_psp_density_valence_fourier(psp, q)
Evaluate the atomic valence charge density in reciprocal space: ρval(q) = ∫{R^3} ρval(r) e^{-iqr} dr = 4π ∫{R+} ρval(r) sin(qr)/qr r^2 dr
DFTK.eval_psp_density_valence_real
— Methodeval_psp_density_valence_real(psp, r)
Evaluate the atomic valence charge density in real space.
DFTK.eval_psp_energy_correction
— Functioneval_psp_energy_correction([T=Float64,] psp, n_electrons)
Evaluate the energy correction to the Ewald electrostatic interaction energy of one unit cell, which is required compared the Ewald expression for point-like nuclei. n_electrons
is the number of electrons per unit cell. This defines the uniform compensating background charge, which is assumed here.
Notice: The returned result is the energy per unit cell and not the energy per volume. To obtain the latter, the caller needs to divide by the unit cell volume.
The energy correction is defined as the limit of the Fourier-transform of the local potential as $q \to 0$, using the same correction as in the Fourier-transform of the local potential: math \lim_{q \to 0} 4π N_{\rm elec} ∫_{ℝ_+} (V(r) - C(r)) \frac{\sin(qr)}{qr} r^2 dr + F[C(r)] = 4π N_{\rm elec} ∫_{ℝ_+} (V(r) + Z/r) r^2 dr
DFTK.eval_psp_local_fourier
— Methodeval_psp_local_fourier(psp, q)
Evaluate the local part of the pseudopotential in reciprocal space:
\[\begin{aligned} V_{\rm loc}(q) &= ∫_{ℝ^3} V_{\rm loc}(r) e^{-iqr} dr \\ &= 4π ∫_{ℝ_+} V_{\rm loc}(r) \frac{\sin(qr)}{q} r dr \end{aligned}\]
In practice, the local potential should be corrected using a Coulomb-like term $C(r) = -Z/r$ to remove the long-range tail of $V_{\rm loc}(r)$ from the integral:
\[\begin{aligned} V_{\rm loc}(q) &= ∫_{ℝ^3} (V_{\rm loc}(r) - C(r)) e^{-iq·r} dr + F[C(r)] \\ &= 4π ∫_{ℝ_+} (V_{\rm loc}(r) + Z/r) \frac{\sin(qr)}{qr} r^2 dr - Z/q^2 -\end{aligned}\]
DFTK.eval_psp_local_real
— Methodeval_psp_local_real(psp, r)
Evaluate the local part of the pseudopotential in real space.
DFTK.eval_psp_projector_fourier
— Methodeval_psp_projector_fourier(psp, i, l, q)
Evaluate the radial part of the i
-th projector for angular momentum l
at the reciprocal vector with modulus q
:
\[\begin{aligned} +\end{aligned}\]
DFTK.eval_psp_local_real
— Methodeval_psp_local_real(psp, r)
Evaluate the local part of the pseudopotential in real space.
DFTK.eval_psp_projector_fourier
— Methodeval_psp_projector_fourier(psp, i, l, q)
Evaluate the radial part of the i
-th projector for angular momentum l
at the reciprocal vector with modulus q
:
\[\begin{aligned} p(q) &= ∫_{ℝ^3} p_{il}(r) e^{-iq·r} dr \\ &= 4π ∫_{ℝ_+} r^2 p_{il}(r) j_l(qr) dr -\end{aligned}\]
DFTK.eval_psp_projector_real
— Methodeval_psp_projector_real(psp, i, l, r)
Evaluate the radial part of the i
-th projector for angular momentum l
in real-space at the vector with modulus r
.
DFTK.filled_occupation
— MethodMaximal occupation of a state (2 for non-spin-polarized electrons, 1 otherwise).
DFTK.find_equivalent_kpt
— MethodFind the equivalent index of the coordinate kcoord
∈ ℝ³ in a list kcoords
∈ [-½, ½)³. ΔG
is the vector of ℤ³ such that kcoords[index] = kcoord + ΔG
.
DFTK.gather_kpts
— MethodGather the distributed data of a quantity depending on k
-Points on the master process and return it. On the other (non-master) processes nothing
is returned.
DFTK.gather_kpts
— MethodGather the distributed $k$-point data on the master process and return it as a PlaneWaveBasis
. On the other (non-master) processes nothing
is returned. The returned object should not be used for computations and only to extract data for post-processing and serialisation to disk.
DFTK.gaussian_valence_charge_density_fourier
— MethodGaussian valence charge density using Abinit's coefficient table, in Fourier space.
DFTK.guess_density
— Functionguess_density(basis::PlaneWaveBasis, method::DensityConstructionMethod,
- magnetic_moments=[]; n_electrons=basis.model.n_electrons)
Build a superposition of atomic densities (SAD) guess density or a rarndom guess density.
The guess atomic densities are taken as one of the following depending on the input method
:
-RandomDensity()
: A random density, normalized to the number of electrons basis.model.n_electrons
. Does not support magnetic moments. -ValenceDensityAuto()
: A combination of the ValenceDensityGaussian
and ValenceDensityPseudo
methods where elements whose pseudopotentials provide numeric valence charge density data use them and elements without use Gaussians. -ValenceDensityGaussian()
: Gaussians of length specified by atom_decay_length
normalized for the correct number of electrons:
\[\hat{ρ}(G) = Z_{\mathrm{valence}} \exp\left(-(2π \text{length} |G|)^2\right)\]
ValenceDensityPseudo()
: Numerical pseudo-atomic valence charge densities from the
pseudopotentials. Will fail if one or more elements in the system has a pseudopotential that does not have valence charge density data.
When magnetic moments are provided, construct a symmetry-broken density guess. The magnetic moments should be specified in units of $μ_B$.
DFTK.hamiltonian_with_total_potential
— MethodReturns a new Hamiltonian with local potential replaced by the given one
DFTK.has_core_density
— MethodCheck presence of model core charge density (non-linear core correction).
DFTK.index_G_vectors
— MethodReturn the index tuple I
such that G_vectors(basis)[I] == G
or the index i
such that G_vectors(basis, kpoint)[i] == G
. Returns nothing if outside the range of valid wave vectors.
DFTK.interpolate_density
— MethodInterpolate a function expressed in a basis basis_in
to a basis basis_out
. This interpolation uses a very basic real-space algorithm, and makes a DWIM-y attempt to take into account the fact that basis_out
can be a supercell of basis_in
.
DFTK.interpolate_kpoint
— MethodInterpolate some data from one $k$-point to another. The interpolation is fast, but not necessarily exact. Intended only to construct guesses for iterative solvers.
DFTK.irfft
— MethodPerform a real valued iFFT; see ifft
. Note that this function silently drops the imaginary part.
DFTK.is_metal
— Methodis_metal(eigenvalues, εF; tol)
Determine whether the provided bands indicate the material is a metal, i.e. where bands are cut by the Fermi level.
DFTK.k_to_kpq_mapping
— MethodReturn the indices of the kpoints
shifted by q
. That is for each kpoint
of the basis
: kpoints[ik].coordinate + q = kpoints[indices[ik]].coordinate
.
DFTK.kcoords_monkhorst_pack
— MethodConstruct the coordinates of the $k$-points in a (shifted) Monkorst-Pack grid
DFTK.kgrid_from_minimal_n_kpoints
— MethodSelects a kgrid size which ensures that at least a n_kpoints
total number of $k$-points are used. The distribution of $k$-points amongst coordinate directions is as uniformly as possible, trying to achieve an identical minimal spacing in all directions.
DFTK.kgrid_from_minimal_spacing
— MethodSelects a kgrid size to ensure a minimal spacing (in inverse Bohrs) between kpoints. A reasonable spacing is 0.13
inverse Bohrs (around $2π * 0.04 \AA^{-1}$).
DFTK.kpath_get_kcoords
— MethodReturn kpoint coordinates in reduced coordinates
DFTK.krange_spin
— MethodReturn the index range of $k$-points that have a particular spin component.
DFTK.list_psp
— Functionlist_psp(element; functional, family, core)
List the pseudopotential files known to DFTK. Allows various ways to restrict the displayed files.
Examples
julia> list_psp(family="hgh")
will list all HGH-type pseudopotentials and
julia> list_psp(family="hgh", functional="lda")
will only list those for LDA (also known as Pade in this context) and
julia> list_psp(:O, core=:semicore)
will list all oxygen semicore pseudopotentials known to DFTK.
DFTK.load_psp
— MethodLoad a pseudopotential file from the library of pseudopotentials. The file is searched in the directory datadir_psp()
and by the key
. If the key
is a path to a valid file, the extension is used to determine the type of the pseudopotential file format and a respective class is returned.
DFTK.load_scfres
— Functionload_scfres(filename)
Load back an scfres
, which has previously been stored with save_scfres
. Note the warning in save_scfres
.
DFTK.model_DFT
— MethodBuild a DFT model from the specified atoms, with the specified functionals.
DFTK.model_LDA
— MethodBuild an LDA model (Perdew & Wang parametrization) from the specified atoms. DOI:10.1103/PhysRevB.45.13244
DFTK.model_PBE
— MethodBuild an PBE-GGA model from the specified atoms. DOI:10.1103/PhysRevLett.77.3865
DFTK.model_SCAN
— MethodBuild a SCAN meta-GGA model from the specified atoms. DOI:10.1103/PhysRevLett.115.036402
DFTK.model_atomic
— MethodConvenience constructor, which builds a standard atomic (kinetic + atomic potential) model. Use extra_terms
to add additional terms.
DFTK.mpi_nprocs
— FunctionNumber of processors used in MPI. Can be called without ensuring initialization.
DFTK.multiply_by_expiqr
— MethodReturn the Fourier coefficients for ψk · e^{i q·r}
in the basis of kpt_out
, where ψk
is defined on a basis kpt_in
.
DFTK.n_elec_core
— MethodReturn the number of core electrons
DFTK.n_elec_valence
— MethodReturn the number of valence electrons
DFTK.n_electrons_from_atoms
— MethodNumber of valence electrons.
DFTK.newton
— Methodnewton(basis::PlaneWaveBasis{T}, ψ0;
+\end{aligned}\]
DFTK.eval_psp_projector_real
— Methodeval_psp_projector_real(psp, i, l, r)
Evaluate the radial part of the i
-th projector for angular momentum l
in real-space at the vector with modulus r
.
DFTK.filled_occupation
— MethodMaximal occupation of a state (2 for non-spin-polarized electrons, 1 otherwise).
DFTK.find_equivalent_kpt
— MethodFind the equivalent index of the coordinate kcoord
∈ ℝ³ in a list kcoords
∈ [-½, ½)³. ΔG
is the vector of ℤ³ such that kcoords[index] = kcoord + ΔG
.
DFTK.gather_kpts
— MethodGather the distributed data of a quantity depending on k
-Points on the master process and return it. On the other (non-master) processes nothing
is returned.
DFTK.gather_kpts
— MethodGather the distributed $k$-point data on the master process and return it as a PlaneWaveBasis
. On the other (non-master) processes nothing
is returned. The returned object should not be used for computations and only to extract data for post-processing and serialisation to disk.
DFTK.gaussian_valence_charge_density_fourier
— MethodGaussian valence charge density using Abinit's coefficient table, in Fourier space.
DFTK.guess_density
— Functionguess_density(basis::PlaneWaveBasis, method::DensityConstructionMethod,
+ magnetic_moments=[]; n_electrons=basis.model.n_electrons)
Build a superposition of atomic densities (SAD) guess density or a rarndom guess density.
The guess atomic densities are taken as one of the following depending on the input method
:
-RandomDensity()
: A random density, normalized to the number of electrons basis.model.n_electrons
. Does not support magnetic moments. -ValenceDensityAuto()
: A combination of the ValenceDensityGaussian
and ValenceDensityPseudo
methods where elements whose pseudopotentials provide numeric valence charge density data use them and elements without use Gaussians. -ValenceDensityGaussian()
: Gaussians of length specified by atom_decay_length
normalized for the correct number of electrons:
\[\hat{ρ}(G) = Z_{\mathrm{valence}} \exp\left(-(2π \text{length} |G|)^2\right)\]
ValenceDensityPseudo()
: Numerical pseudo-atomic valence charge densities from the
pseudopotentials. Will fail if one or more elements in the system has a pseudopotential that does not have valence charge density data.
When magnetic moments are provided, construct a symmetry-broken density guess. The magnetic moments should be specified in units of $μ_B$.
DFTK.hamiltonian_with_total_potential
— MethodReturns a new Hamiltonian with local potential replaced by the given one
DFTK.has_core_density
— MethodCheck presence of model core charge density (non-linear core correction).
DFTK.index_G_vectors
— MethodReturn the index tuple I
such that G_vectors(basis)[I] == G
or the index i
such that G_vectors(basis, kpoint)[i] == G
. Returns nothing if outside the range of valid wave vectors.
DFTK.interpolate_density
— MethodInterpolate a function expressed in a basis basis_in
to a basis basis_out
. This interpolation uses a very basic real-space algorithm, and makes a DWIM-y attempt to take into account the fact that basis_out
can be a supercell of basis_in
.
DFTK.interpolate_kpoint
— MethodInterpolate some data from one $k$-point to another. The interpolation is fast, but not necessarily exact. Intended only to construct guesses for iterative solvers.
DFTK.irfft
— MethodPerform a real valued iFFT; see ifft
. Note that this function silently drops the imaginary part.
DFTK.is_metal
— Methodis_metal(eigenvalues, εF; tol)
Determine whether the provided bands indicate the material is a metal, i.e. where bands are cut by the Fermi level.
DFTK.k_to_kpq_mapping
— MethodReturn the indices of the kpoints
shifted by q
. That is for each kpoint
of the basis
: kpoints[ik].coordinate + q = kpoints[indices[ik]].coordinate
.
DFTK.kcoords_monkhorst_pack
— MethodConstruct the coordinates of the $k$-points in a (shifted) Monkorst-Pack grid
DFTK.kgrid_from_minimal_n_kpoints
— MethodSelects a kgrid size which ensures that at least a n_kpoints
total number of $k$-points are used. The distribution of $k$-points amongst coordinate directions is as uniformly as possible, trying to achieve an identical minimal spacing in all directions.
DFTK.kgrid_from_minimal_spacing
— MethodSelects a kgrid size to ensure a minimal spacing (in inverse Bohrs) between kpoints. A reasonable spacing is 0.13
inverse Bohrs (around $2π * 0.04 \AA^{-1}$).
DFTK.kpath_get_kcoords
— MethodReturn kpoint coordinates in reduced coordinates
DFTK.krange_spin
— MethodReturn the index range of $k$-points that have a particular spin component.
DFTK.list_psp
— Functionlist_psp(element; functional, family, core)
List the pseudopotential files known to DFTK. Allows various ways to restrict the displayed files.
Examples
julia> list_psp(family="hgh")
will list all HGH-type pseudopotentials and
julia> list_psp(family="hgh", functional="lda")
will only list those for LDA (also known as Pade in this context) and
julia> list_psp(:O, core=:semicore)
will list all oxygen semicore pseudopotentials known to DFTK.
DFTK.load_psp
— MethodLoad a pseudopotential file from the library of pseudopotentials. The file is searched in the directory datadir_psp()
and by the key
. If the key
is a path to a valid file, the extension is used to determine the type of the pseudopotential file format and a respective class is returned.
DFTK.load_scfres
— Functionload_scfres(filename)
Load back an scfres
, which has previously been stored with save_scfres
. Note the warning in save_scfres
.
DFTK.model_DFT
— MethodBuild a DFT model from the specified atoms, with the specified functionals.
DFTK.model_LDA
— MethodBuild an LDA model (Perdew & Wang parametrization) from the specified atoms. DOI:10.1103/PhysRevB.45.13244
DFTK.model_PBE
— MethodBuild an PBE-GGA model from the specified atoms. DOI:10.1103/PhysRevLett.77.3865
DFTK.model_SCAN
— MethodBuild a SCAN meta-GGA model from the specified atoms. DOI:10.1103/PhysRevLett.115.036402
DFTK.model_atomic
— MethodConvenience constructor, which builds a standard atomic (kinetic + atomic potential) model. Use extra_terms
to add additional terms.
DFTK.mpi_nprocs
— FunctionNumber of processors used in MPI. Can be called without ensuring initialization.
DFTK.multiply_by_expiqr
— MethodReturn the Fourier coefficients for ψk · e^{i q·r}
in the basis of kpt_out
, where ψk
is defined on a basis kpt_in
.
DFTK.n_elec_core
— MethodReturn the number of core electrons
DFTK.n_elec_valence
— MethodReturn the number of valence electrons
DFTK.n_electrons_from_atoms
— MethodNumber of valence electrons.
DFTK.newton
— Methodnewton(basis::PlaneWaveBasis{T}, ψ0;
tol=1e-6, tol_cg=tol / 100, maxiter=20, callback=ScfDefaultCallback(),
- is_converged=ScfConvergenceDensity(tol))
Newton algorithm. Be careful that the starting point needs to be not too far from the solution.
DFTK.next_compatible_fft_size
— MethodFind the next compatible FFT size Sizes must (a) be a product of small primes only and (b) contain the factors. If smallprimes is empty (a) is skipped.
DFTK.next_density
— FunctionObtain new density ρ by diagonalizing ham
. Follows the policy imposed by the bands
data structure to determine and adjust the number of bands to be computed.
DFTK.norm2
— MethodSquare of the ℓ²-norm.
DFTK.norm_cplx
— MethodComplex-analytic extension of LinearAlgebra.norm(x)
from real to complex inputs. Needed for phonons as we want to perform a matrix-vector product f'(x)·h
, where f
is a real-to-real function and h
a complex vector. To do this using automatic differentiation, we can extend analytically f to accept complex inputs, then differentiate t -> f(x+t·h)
. This will fail if non-analytic functions like norm are used for complex inputs, and therefore we have to redefine it.
DFTK.normalize_kpoint_coordinate
— MethodBring $k$-point coordinates into the range [-0.5, 0.5)
DFTK.overlap_Mmn_k_kpb
— MethodComputes the matrix $[M^{k,b}]_{m,n} = \langle u_{m,k} | u_{n,k+b} \rangle$ for given k
, kpb
= $k+b$.
G_shift
is the "shifting" vector, correction due to the periodicity conditions imposed on $k \to ψ_k$. It is non zero if kpb
is taken in another unit cell of the reciprocal lattice. We use here that: $u_{n(k + G_{\rm shift})}(r) = e^{-i*\langle G_{\rm shift},r \rangle} u_{nk}$.
DFTK.plot_bandstructure
— FunctionCompute and plot the band structure. n_bands
selects the number of bands to compute. If this value is absent and an scfres
is used to start the calculation a default of n_bands_scf + 5sqrt(n_bands_scf)
is used. The unit used to plot the bands can be selected using the unit
parameter. Like in the rest of DFTK Hartree is used by default. Another standard choices is unit=u"eV"
(electron volts). The kline_density
is given in number of $k$-points per inverse bohrs (i.e. overall in units of length).
DFTK.plot_dos
— FunctionPlot the density of states over a reasonable range. Requires to load Plots.jl
beforehand.
DFTK.psp_local_polynomial
— FunctionThe local potential of a HGH pseudopotentials in reciprocal space can be brought to the form $Q(t) / (t^2 exp(t^2 / 2))$ where $t = r_\text{loc} q$ and Q
is a polynomial of at most degree 8. This function returns Q
.
DFTK.psp_projector_polynomial
— FunctionThe nonlocal projectors of a HGH pseudopotentials in reciprocal space can be brought to the form $Q(t) exp(-t^2 / 2)$ where $t = r_l q$ and Q
is a polynomial. This function returns Q
.
DFTK.qcut_psp_local
— MethodEstimate an upper bound for the argument q
after which abs(eval_psp_local_fourier(psp, q))
is a strictly decreasing function.
DFTK.qcut_psp_projector
— MethodEstimate an upper bound for the argument q
after which eval_psp_projector_fourier(psp, q)
is a strictly decreasing function.
DFTK.r_vectors
— Methodr_vectors(basis::PlaneWaveBasis)
The list of $r$ vectors, in reduced coordinates. By convention, this is in [0,1)^3.
DFTK.r_vectors_cart
— Methodr_vectors_cart(basis::PlaneWaveBasis)
The list of $r$ vectors, in cartesian coordinates.
DFTK.random_density
— MethodBuild a random charge density normalized to the provided number of electrons.
DFTK.read_w90_nnkp
— MethodRead the .nnkp file provided by the preprocessing routine of Wannier90 (i.e. "wannier90.x -pp prefix") Returns:
- the array 'nnkpts' of k points, their respective nearest neighbors and associated shifing vectors (non zero if the neighbor is located in another cell).
- the number 'nntot' of neighbors per k point.
TODO: add the possibility to exclude bands
DFTK.run_wannier90
— FunctionWannerize the obtained bands using wannier90. By default all converged bands from the scfres
are employed (change with n_bands
kwargs) and n_wannier = n_bands
wannier functions are computed using random Gaussians as guesses. All keyword arguments supported by Wannier90 for the disentanglement may be added as keyword arguments. The function returns the fileprefix
.
Currently this is an experimental feature, which has not yet been tested to full depth. The interface is considered unstable and may change incompatibly in the future. Use at your own risk and please report bugs in case you encounter any.
DFTK.save_scfres
— Methodsave_scfres(filename, scfres)
Save an scfres
obtained from self_consistent_field
to a file. The format is determined from the file extension. Currently the following file extensions are recognized and supported:
- jld2: A JLD2 file. Stores the complete state and can be used (with
load_scfres
) to restart an SCF from a checkpoint or post-process an SCF solution. See Saving SCF results on disk and SCF checkpoints for details. - vts: A VTK file for visualisation e.g. in paraview. Stores the density, spin density and some metadata (energy, Fermi level, occupation etc.). Supports these keyword arguments:
save_ψ
: Save the real-space representation of the orbitals as well (may lead to larger files).extra_data
:Dict{String,Array}
with additional data on the 3D real-space grid to store into the VTK file.
- json: A JSON file with basic information about the SCF run. Stores for example the number of iterations, occupations, norm of the most recent density change, eigenvalues, Fermi level etc.
No guarantees are made with respect to this function at this point. It may change incompatibly between DFTK versions or stop working / be removed in the future.
DFTK.scf_anderson_solver
— FunctionCreate a simple anderson-accelerated SCF solver. m
specifies the number of steps to keep the history of.
DFTK.scf_damping_quadratic_model
— MethodUse the two iteration states info
and info_next
to find a damping value from a quadratic model for the SCF energy. Returns nothing
if the constructed model is not considered trustworthy, else returns the suggested damping.
DFTK.scf_damping_solver
— FunctionCreate a damped SCF solver updating the density as x = β * x_new + (1 - β) * x
DFTK.select_eigenpairs_all_kblocks
— MethodFunction to select a subset of eigenpairs on each $k$-Point. Works on the Tuple returned by diagonalize_all_kblocks
.
DFTK.self_consistent_field
— Methodself_consistent_field(basis; [tol, mixing, damping, ρ, ψ])
Solve the Kohn-Sham equations with a density-based SCF algorithm using damped, preconditioned iterations where $ρ_\text{next} = α P^{-1} (ρ_\text{out} - ρ_\text{in})$.
Overview of parameters:
ρ
: Initial densityψ
: Initial orbitalstol
: Tolerance for the density change ($\|ρ_\text{out} - ρ_\text{in}\|$) to flag convergence. Default is1e-6
.is_converged
: Convergence control callback. Typical objects passed here areDFTK.ScfConvergenceDensity(tol)
(the default),DFTK.ScfConvergenceEnergy(tol)
orDFTK.ScfConvergenceForce(tol)
.maxiter
: Maximal number of SCF iterationsmixing
: Mixing method, which determines the preconditioner $P^{-1}$ in the above equation. Typical mixings areLdosMixing
,KerkerMixing
,SimpleMixing
orDielectricMixing
. Default isLdosMixing()
damping
: Damping parameter $α$ in the above equation. Default is0.8
.nbandsalg
: By default DFTK usesnbandsalg=AdaptiveBands(model)
, which adaptively determines the number of bands to compute. If you want to influence this algorithm or use a predefined number of bands in each SCF step, pass aFixedBands
orAdaptiveBands
. Beware that with non-zero temperature, the convergence of the SCF algorithm may be limited by thedefault_occupation_threshold
parameter. For highly accurate calculations we thus recommend increasing thedefault_occupation_threshold
of theAdaptiveBands
.callback
: Function called at each SCF iteration. Usually takes care of printing the intermediate state.
DFTK.sin2pi
— MethodFunction to compute sin(2π x)
DFTK.solve_ΩplusK
— Methodsolve_ΩplusK(basis::PlaneWaveBasis{T}, ψ, res, occupation;
- tol=1e-10, verbose=false) where {T}
Return δψ where (Ω+K) δψ = rhs
DFTK.solve_ΩplusK_split
— MethodSolve the problem (Ω+K) δψ = rhs
using a split algorithm, where rhs
is typically -δHextψ
(the negative matvec of an external perturbation with the SCF orbitals ψ
) and δψ
is the corresponding total variation in the orbitals ψ
. Additionally returns: - δρ
: Total variation in density) - δHψ
: Total variation in Hamiltonian applied to orbitals - δeigenvalues
: Total variation in eigenvalues - δVind
: Change in potential induced by δρ
(the term needed on top of δHextψ
to get δHψ
).
DFTK.spglib_standardize_cell
— MethodReturns crystallographic conventional cell according to the International Table of Crystallography Vol A (ITA) in case primitive=false
. If primitive=true
the primitive lattice is returned in the convention of the reference work of Cracknell, Davies, Miller, and Love (CDML). Of note this has minor differences to the primitive setting choice made in the ITA.
DFTK.sphericalbesselj_fast
— Methodsphericalbesselj_fast(l::Integer, x::Number)
Returns the spherical Bessel function of the first kind jl(x). Consistent with https://en.wikipedia.org/wiki/Besselfunction#SphericalBesselfunctions and with SpecialFunctions.sphericalbesselj
. Specialized for integer 0 <= l <= 5
.
DFTK.spin_components
— MethodExplicit spin components of the KS orbitals and the density
DFTK.split_evenly
— MethodSplit an iterable evenly into N chunks, which will be returned.
DFTK.standardize_atoms
— FunctionApply various standardisations to a lattice and a list of atoms. It uses spglib to detect symmetries (within tol_symmetry
), then cleans up the lattice according to the symmetries (unless correct_symmetry
is false
) and returns the resulting standard lattice and atoms. If primitive
is true
(default) the primitive unit cell is returned, else the conventional unit cell is returned.
DFTK.symmetries_preserving_kgrid
— MethodFilter out the symmetry operations that don't respect the symmetries of the discrete BZ grid
DFTK.symmetries_preserving_rgrid
— MethodFilter out the symmetry operations that don't respect the symmetries of the discrete real-space grid
DFTK.symmetrize_forces
— MethodSymmetrize the forces in reduced coordinates, forces given as an array forces[iel][α,i]
DFTK.symmetrize_stresses
— MethodSymmetrize the stress tensor, given as a Matrix in cartesian coordinates
DFTK.symmetrize_ρ
— MethodSymmetrize a density by applying all the basis (by default) symmetries and forming the average.
DFTK.symmetry_operations
— FunctionReturn the symmetries given an atomic structure with optionally designated magnetic moments on each of the atoms. The symmetries are determined using spglib.
DFTK.symmetry_operations
— MethodReturn the Symmetry operations given a hall_number
.
This function allows to directly access to the space group operations in the spglib
database. To specify the space group type with a specific choice, hall_number
is used.
The definition of hall_number
is found at Space group type.
DFTK.synchronize_device
— MethodSynchronize data and finish all operations on the execution stream of the device. This needs to be called explicitly before a task finishes (e.g. in an @spawn
block).
DFTK.to_cpu
— MethodTransfer an array from a device (typically a GPU) to the CPU.
DFTK.to_device
— MethodTransfer an array to a particular device (typically a GPU)
DFTK.todict
— MethodConvert an Energies
struct to a dictionary representation
DFTK.total_local_potential
— MethodGet the total local potential of the given Hamiltonian, in real space in the spin components.
DFTK.transfer_blochwave
— MethodTransfer Bloch wave between two basis sets. Limited feature set.
DFTK.transfer_blochwave_kpt
— MethodTransfer an array ψk_in
expanded on kpt_in
, and produce $ψ(r) e^{i ΔG·r}$ expanded on kpt_out
. It is mostly useful for phonons. Beware: ψk_out
can lose information if the shift ΔG
is large or if the G_vectors
differ between k
-points.
DFTK.transfer_blochwave_kpt
— MethodTransfer an array ψk
defined on basisin $k$-point kptin to basisout $k$-point kptout.
DFTK.transfer_density
— MethodTransfer density (in real space) between two basis sets.
This function is fast by transferring only the Fourier coefficients from the small basis to the big basis.
Note that this implies that for even-sized small FFT grids doing the transfer small -> big -> small is not an identity (as the small basis has an unmatched Fourier component and the identity $c_G = c_{-G}^\ast$ does not fully hold).
Note further that for the direction big -> small employing this function does not give the same answer as using first transfer_blochwave
and then compute_density
.
DFTK.transfer_mapping
— MethodCompute the index mapping between the global grids of two bases. Returns an iterator of 8 pairs (block_in, block_out)
. Iterated over these pairs x_out_fourier[block_out, :] = x_in_fourier[block_in, :]
does the transfer from the Fourier coefficients x_in_fourier
(defined on basis_in
) to x_out_fourier
(defined on basis_out
, equally provided as Fourier coefficients).
DFTK.transfer_mapping
— MethodCompute the index mapping between two bases. Returns two arrays idcs_in
and idcs_out
such that ψkout[idcs_out] = ψkin[idcs_in]
does the transfer from ψkin
(defined on basis_in
and kpt_in
) to ψkout
(defined on basis_out
and kpt_out
).
DFTK.unfold_bz
— Method" Convert a basis
into one that doesn't use BZ symmetry. This is mainly useful for debug purposes (e.g. in cases we don't want to bother thinking about symmetries).
DFTK.versioninfo
— FunctionDFTK.versioninfo([io::IO=stdout])
Summary of version and configuration of DFTK and its key dependencies.
DFTK.weighted_ksum
— MethodSum an array over kpoints, taking weights into account
DFTK.write_w90_eig
— MethodWrite the eigenvalues in a format readable by Wannier90.
DFTK.write_w90_win
— MethodWrite a win file at the indicated prefix. Parameters to Wannier90 can be added as kwargs: e.g. num_iter=500
.
DFTK.ylm_real
— MethodReturns the (l,m) real spherical harmonic Ylm(r). Consistent with https://en.wikipedia.org/wiki/Tableofsphericalharmonics#Realsphericalharmonics
DFTK.zeros_like
— FunctionCreate an array of same "array type" as X filled with zeros, minimizing the number of allocations. This unifies CPU and GPU code, as the output will always be on the same device as the input.
DFTK.@timing
— MacroShortened version of the @timeit
macro from TimerOutputs
, which writes to the DFTK timer.
DFTK.Smearing.A
— MethodA
term in the Hermite delta expansion
DFTK.Smearing.H
— MethodStandard Hermite function using physicist's convention.
DFTK.Smearing.entropy
— MethodEntropy. Note that this is a function of the energy x
, not of occupation(x)
. This function satisfies s' = x f' (see https://www.vasp.at/vasp-workshop/k-points.pdf p. 12 and https://arxiv.org/pdf/1805.07144.pdf p. 18.
DFTK.Smearing.occupation
— Functionoccupation(S::SmearingFunction, x)
Occupation at x
, where in practice x = (ε - εF) / temperature
. If temperature is zero, (ε-εF)/temperature = ±∞
. The occupation function is required to give 1 and 0 respectively in these cases.
DFTK.Smearing.occupation_derivative
— MethodDerivative of the occupation function, approximation to minus the delta function.
DFTK.Smearing.occupation_divided_difference
— Method(f(x) - f(y))/(x - y), computed stably in the case where x and y are close
Settings
This document was generated with Documenter.jl version 1.1.2 on Tuesday 7 November 2023. Using Julia version 1.9.3.
Newton algorithm. Be careful that the starting point needs to be not too far from the solution.
DFTK.next_compatible_fft_size
— MethodFind the next compatible FFT size Sizes must (a) be a product of small primes only and (b) contain the factors. If smallprimes is empty (a) is skipped.
DFTK.next_density
— FunctionObtain new density ρ by diagonalizing ham
. Follows the policy imposed by the bands
data structure to determine and adjust the number of bands to be computed.
DFTK.norm2
— MethodSquare of the ℓ²-norm.
DFTK.norm_cplx
— MethodComplex-analytic extension of LinearAlgebra.norm(x)
from real to complex inputs. Needed for phonons as we want to perform a matrix-vector product f'(x)·h
, where f
is a real-to-real function and h
a complex vector. To do this using automatic differentiation, we can extend analytically f to accept complex inputs, then differentiate t -> f(x+t·h)
. This will fail if non-analytic functions like norm are used for complex inputs, and therefore we have to redefine it.
DFTK.normalize_kpoint_coordinate
— MethodBring $k$-point coordinates into the range [-0.5, 0.5)
DFTK.overlap_Mmn_k_kpb
— MethodComputes the matrix $[M^{k,b}]_{m,n} = \langle u_{m,k} | u_{n,k+b} \rangle$ for given k
, kpb
= $k+b$.
G_shift
is the "shifting" vector, correction due to the periodicity conditions imposed on $k \to ψ_k$. It is non zero if kpb
is taken in another unit cell of the reciprocal lattice. We use here that: $u_{n(k + G_{\rm shift})}(r) = e^{-i*\langle G_{\rm shift},r \rangle} u_{nk}$.
DFTK.plot_bandstructure
— FunctionCompute and plot the band structure. n_bands
selects the number of bands to compute. If this value is absent and an scfres
is used to start the calculation a default of n_bands_scf + 5sqrt(n_bands_scf)
is used. The unit used to plot the bands can be selected using the unit
parameter. Like in the rest of DFTK Hartree is used by default. Another standard choices is unit=u"eV"
(electron volts). The kline_density
is given in number of $k$-points per inverse bohrs (i.e. overall in units of length).
DFTK.plot_dos
— FunctionPlot the density of states over a reasonable range. Requires to load Plots.jl
beforehand.
DFTK.psp_local_polynomial
— FunctionThe local potential of a HGH pseudopotentials in reciprocal space can be brought to the form $Q(t) / (t^2 exp(t^2 / 2))$ where $t = r_\text{loc} q$ and Q
is a polynomial of at most degree 8. This function returns Q
.
DFTK.psp_projector_polynomial
— FunctionThe nonlocal projectors of a HGH pseudopotentials in reciprocal space can be brought to the form $Q(t) exp(-t^2 / 2)$ where $t = r_l q$ and Q
is a polynomial. This function returns Q
.
DFTK.qcut_psp_local
— MethodEstimate an upper bound for the argument q
after which abs(eval_psp_local_fourier(psp, q))
is a strictly decreasing function.
DFTK.qcut_psp_projector
— MethodEstimate an upper bound for the argument q
after which eval_psp_projector_fourier(psp, q)
is a strictly decreasing function.
DFTK.r_vectors
— Methodr_vectors(basis::PlaneWaveBasis)
The list of $r$ vectors, in reduced coordinates. By convention, this is in [0,1)^3.
DFTK.r_vectors_cart
— Methodr_vectors_cart(basis::PlaneWaveBasis)
The list of $r$ vectors, in cartesian coordinates.
DFTK.random_density
— MethodBuild a random charge density normalized to the provided number of electrons.
DFTK.read_w90_nnkp
— MethodRead the .nnkp file provided by the preprocessing routine of Wannier90 (i.e. "wannier90.x -pp prefix") Returns:
- the array 'nnkpts' of k points, their respective nearest neighbors and associated shifing vectors (non zero if the neighbor is located in another cell).
- the number 'nntot' of neighbors per k point.
TODO: add the possibility to exclude bands
DFTK.run_wannier90
— FunctionWannerize the obtained bands using wannier90. By default all converged bands from the scfres
are employed (change with n_bands
kwargs) and n_wannier = n_bands
wannier functions are computed using random Gaussians as guesses. All keyword arguments supported by Wannier90 for the disentanglement may be added as keyword arguments. The function returns the fileprefix
.
Currently this is an experimental feature, which has not yet been tested to full depth. The interface is considered unstable and may change incompatibly in the future. Use at your own risk and please report bugs in case you encounter any.
DFTK.save_scfres
— Methodsave_scfres(filename, scfres)
Save an scfres
obtained from self_consistent_field
to a file. The format is determined from the file extension. Currently the following file extensions are recognized and supported:
- jld2: A JLD2 file. Stores the complete state and can be used (with
load_scfres
) to restart an SCF from a checkpoint or post-process an SCF solution. See Saving SCF results on disk and SCF checkpoints for details. - vts: A VTK file for visualisation e.g. in paraview. Stores the density, spin density and some metadata (energy, Fermi level, occupation etc.). Supports these keyword arguments:
save_ψ
: Save the real-space representation of the orbitals as well (may lead to larger files).extra_data
:Dict{String,Array}
with additional data on the 3D real-space grid to store into the VTK file.
- json: A JSON file with basic information about the SCF run. Stores for example the number of iterations, occupations, norm of the most recent density change, eigenvalues, Fermi level etc.
No guarantees are made with respect to this function at this point. It may change incompatibly between DFTK versions or stop working / be removed in the future.
DFTK.scf_anderson_solver
— FunctionCreate a simple anderson-accelerated SCF solver. m
specifies the number of steps to keep the history of.
DFTK.scf_damping_quadratic_model
— MethodUse the two iteration states info
and info_next
to find a damping value from a quadratic model for the SCF energy. Returns nothing
if the constructed model is not considered trustworthy, else returns the suggested damping.
DFTK.scf_damping_solver
— FunctionCreate a damped SCF solver updating the density as x = β * x_new + (1 - β) * x
DFTK.select_eigenpairs_all_kblocks
— MethodFunction to select a subset of eigenpairs on each $k$-Point. Works on the Tuple returned by diagonalize_all_kblocks
.
DFTK.self_consistent_field
— Methodself_consistent_field(basis; [tol, mixing, damping, ρ, ψ])
Solve the Kohn-Sham equations with a density-based SCF algorithm using damped, preconditioned iterations where $ρ_\text{next} = α P^{-1} (ρ_\text{out} - ρ_\text{in})$.
Overview of parameters:
ρ
: Initial densityψ
: Initial orbitalstol
: Tolerance for the density change ($\|ρ_\text{out} - ρ_\text{in}\|$) to flag convergence. Default is1e-6
.is_converged
: Convergence control callback. Typical objects passed here areDFTK.ScfConvergenceDensity(tol)
(the default),DFTK.ScfConvergenceEnergy(tol)
orDFTK.ScfConvergenceForce(tol)
.maxiter
: Maximal number of SCF iterationsmixing
: Mixing method, which determines the preconditioner $P^{-1}$ in the above equation. Typical mixings areLdosMixing
,KerkerMixing
,SimpleMixing
orDielectricMixing
. Default isLdosMixing()
damping
: Damping parameter $α$ in the above equation. Default is0.8
.nbandsalg
: By default DFTK usesnbandsalg=AdaptiveBands(model)
, which adaptively determines the number of bands to compute. If you want to influence this algorithm or use a predefined number of bands in each SCF step, pass aFixedBands
orAdaptiveBands
. Beware that with non-zero temperature, the convergence of the SCF algorithm may be limited by thedefault_occupation_threshold
parameter. For highly accurate calculations we thus recommend increasing thedefault_occupation_threshold
of theAdaptiveBands
.callback
: Function called at each SCF iteration. Usually takes care of printing the intermediate state.
DFTK.sin2pi
— MethodFunction to compute sin(2π x)
DFTK.solve_ΩplusK
— Methodsolve_ΩplusK(basis::PlaneWaveBasis{T}, ψ, res, occupation;
+ tol=1e-10, verbose=false) where {T}
Return δψ where (Ω+K) δψ = rhs
DFTK.solve_ΩplusK_split
— MethodSolve the problem (Ω+K) δψ = rhs
using a split algorithm, where rhs
is typically -δHextψ
(the negative matvec of an external perturbation with the SCF orbitals ψ
) and δψ
is the corresponding total variation in the orbitals ψ
. Additionally returns: - δρ
: Total variation in density) - δHψ
: Total variation in Hamiltonian applied to orbitals - δeigenvalues
: Total variation in eigenvalues - δVind
: Change in potential induced by δρ
(the term needed on top of δHextψ
to get δHψ
).
DFTK.spglib_standardize_cell
— MethodReturns crystallographic conventional cell according to the International Table of Crystallography Vol A (ITA) in case primitive=false
. If primitive=true
the primitive lattice is returned in the convention of the reference work of Cracknell, Davies, Miller, and Love (CDML). Of note this has minor differences to the primitive setting choice made in the ITA.
DFTK.sphericalbesselj_fast
— Methodsphericalbesselj_fast(l::Integer, x::Number)
Returns the spherical Bessel function of the first kind jl(x). Consistent with https://en.wikipedia.org/wiki/Besselfunction#SphericalBesselfunctions and with SpecialFunctions.sphericalbesselj
. Specialized for integer 0 <= l <= 5
.
DFTK.spin_components
— MethodExplicit spin components of the KS orbitals and the density
DFTK.split_evenly
— MethodSplit an iterable evenly into N chunks, which will be returned.
DFTK.standardize_atoms
— FunctionApply various standardisations to a lattice and a list of atoms. It uses spglib to detect symmetries (within tol_symmetry
), then cleans up the lattice according to the symmetries (unless correct_symmetry
is false
) and returns the resulting standard lattice and atoms. If primitive
is true
(default) the primitive unit cell is returned, else the conventional unit cell is returned.
DFTK.symmetries_preserving_kgrid
— MethodFilter out the symmetry operations that don't respect the symmetries of the discrete BZ grid
DFTK.symmetries_preserving_rgrid
— MethodFilter out the symmetry operations that don't respect the symmetries of the discrete real-space grid
DFTK.symmetrize_forces
— MethodSymmetrize the forces in reduced coordinates, forces given as an array forces[iel][α,i]
DFTK.symmetrize_stresses
— MethodSymmetrize the stress tensor, given as a Matrix in cartesian coordinates
DFTK.symmetrize_ρ
— MethodSymmetrize a density by applying all the basis (by default) symmetries and forming the average.
DFTK.symmetry_operations
— FunctionReturn the symmetries given an atomic structure with optionally designated magnetic moments on each of the atoms. The symmetries are determined using spglib.
DFTK.symmetry_operations
— MethodReturn the Symmetry operations given a hall_number
.
This function allows to directly access to the space group operations in the spglib
database. To specify the space group type with a specific choice, hall_number
is used.
The definition of hall_number
is found at Space group type.
DFTK.synchronize_device
— MethodSynchronize data and finish all operations on the execution stream of the device. This needs to be called explicitly before a task finishes (e.g. in an @spawn
block).
DFTK.to_cpu
— MethodTransfer an array from a device (typically a GPU) to the CPU.
DFTK.to_device
— MethodTransfer an array to a particular device (typically a GPU)
DFTK.todict
— MethodConvert an Energies
struct to a dictionary representation
DFTK.total_local_potential
— MethodGet the total local potential of the given Hamiltonian, in real space in the spin components.
DFTK.transfer_blochwave
— MethodTransfer Bloch wave between two basis sets. Limited feature set.
DFTK.transfer_blochwave_kpt
— MethodTransfer an array ψk_in
expanded on kpt_in
, and produce $ψ(r) e^{i ΔG·r}$ expanded on kpt_out
. It is mostly useful for phonons. Beware: ψk_out
can lose information if the shift ΔG
is large or if the G_vectors
differ between k
-points.
DFTK.transfer_blochwave_kpt
— MethodTransfer an array ψk
defined on basisin $k$-point kptin to basisout $k$-point kptout.
DFTK.transfer_density
— MethodTransfer density (in real space) between two basis sets.
This function is fast by transferring only the Fourier coefficients from the small basis to the big basis.
Note that this implies that for even-sized small FFT grids doing the transfer small -> big -> small is not an identity (as the small basis has an unmatched Fourier component and the identity $c_G = c_{-G}^\ast$ does not fully hold).
Note further that for the direction big -> small employing this function does not give the same answer as using first transfer_blochwave
and then compute_density
.
DFTK.transfer_mapping
— MethodCompute the index mapping between the global grids of two bases. Returns an iterator of 8 pairs (block_in, block_out)
. Iterated over these pairs x_out_fourier[block_out, :] = x_in_fourier[block_in, :]
does the transfer from the Fourier coefficients x_in_fourier
(defined on basis_in
) to x_out_fourier
(defined on basis_out
, equally provided as Fourier coefficients).
DFTK.transfer_mapping
— MethodCompute the index mapping between two bases. Returns two arrays idcs_in
and idcs_out
such that ψkout[idcs_out] = ψkin[idcs_in]
does the transfer from ψkin
(defined on basis_in
and kpt_in
) to ψkout
(defined on basis_out
and kpt_out
).
DFTK.unfold_bz
— Method" Convert a basis
into one that doesn't use BZ symmetry. This is mainly useful for debug purposes (e.g. in cases we don't want to bother thinking about symmetries).
DFTK.versioninfo
— FunctionDFTK.versioninfo([io::IO=stdout])
Summary of version and configuration of DFTK and its key dependencies.
DFTK.weighted_ksum
— MethodSum an array over kpoints, taking weights into account
DFTK.write_w90_eig
— MethodWrite the eigenvalues in a format readable by Wannier90.
DFTK.write_w90_win
— MethodWrite a win file at the indicated prefix. Parameters to Wannier90 can be added as kwargs: e.g. num_iter=500
.
DFTK.ylm_real
— MethodReturns the (l,m) real spherical harmonic Ylm(r). Consistent with https://en.wikipedia.org/wiki/Tableofsphericalharmonics#Realsphericalharmonics
DFTK.zeros_like
— FunctionCreate an array of same "array type" as X filled with zeros, minimizing the number of allocations. This unifies CPU and GPU code, as the output will always be on the same device as the input.
DFTK.@timing
— MacroShortened version of the @timeit
macro from TimerOutputs
, which writes to the DFTK timer.
DFTK.Smearing.A
— MethodA
term in the Hermite delta expansion
DFTK.Smearing.H
— MethodStandard Hermite function using physicist's convention.
DFTK.Smearing.entropy
— MethodEntropy. Note that this is a function of the energy x
, not of occupation(x)
. This function satisfies s' = x f' (see https://www.vasp.at/vasp-workshop/k-points.pdf p. 12 and https://arxiv.org/pdf/1805.07144.pdf p. 18.
DFTK.Smearing.occupation
— Functionoccupation(S::SmearingFunction, x)
Occupation at x
, where in practice x = (ε - εF) / temperature
. If temperature is zero, (ε-εF)/temperature = ±∞
. The occupation function is required to give 1 and 0 respectively in these cases.
DFTK.Smearing.occupation_derivative
— MethodDerivative of the occupation function, approximation to minus the delta function.
DFTK.Smearing.occupation_divided_difference
— Method(f(x) - f(y))/(x - y), computed stably in the case where x and y are close