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.AdaptiveBandsType

Dynamically 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.

source
DFTK.AdaptiveDiagtolType

Algorithm for the tolerance used for the next diagonalization. This function takes $|ρ_{\rm next} - ρ_{\rm 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.

source
DFTK.AtomicNonlocalType

Nonlocal term coming from norm-conserving pseudopotentials in Kleinmann-Bylander form.

\[\text{Energy} = ∑_a ∑_{ij} ∑_{n} f_n \braket{ψ_n}{{\rm proj}_{ai}} D_{ij} \braket{{\rm proj}_{aj}}{ψ_n}.\]

source
DFTK.BlowupCHVType

Blow-up function as proposed in arXiv: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.

source
DFTK.DFTKCalculatorMethod
DFTKCalculator(
    params::DFTK.DFTKParameters
) -> DFTKCalculator

Construct a AtomsCalculators compatible calculator for DFTK. The model_kwargs are passed onto the Model constructor, the basis_kwargs to the PlaneWaveBasis constructor, the scf_kwargs to self_consistent_field. At the very least the DFT functionals and the Ecut needs to be specified.

By default the calculator preserves the symmetries that are stored inside the state (the basis is re-built, but symmetries are fixed and not re-computed).

Example

julia> DFTKCalculator(; model_kwargs=(; functionals=[:lda_x, :lda_c_vwn]),
                        basis_kwargs=(; Ecut=10, kgrid=(2, 2, 2)),
                        scf_kwargs=(; tol=1e-4))
source
DFTK.DielectricMixingType

We use a simplification of the Resta model 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.

source
DFTK.DielectricModelType

A 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.

source
DFTK.DivAgradOperatorType

Nonlocal "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$).

source
DFTK.ElementCohenBergstresserMethod
ElementCohenBergstresser(
    key;
    lattice_constant
) -> ElementCohenBergstresser

Element 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")

source
DFTK.ElementCoulombMethod
ElementCoulomb(key; mass) -> ElementCoulomb

Element 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")

source
DFTK.ElementGaussianMethod
ElementGaussian(α, L; symbol) -> ElementGaussian

Element interacting with electrons via a Gaussian potential. Symbol is non-mandatory.

source
DFTK.ElementPspMethod
ElementPsp(key; psp, mass)

Element 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")

source
DFTK.EnergiesType

A simple struct to contain a vector of energies, and utilities to print them in a nice format.

source
DFTK.EntropyType

Entropy 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.

source
DFTK.EwaldType

Ewald 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.

source
DFTK.ExplicitKpointsType

Explicitly define the k-points along which to perform BZ sampling. (Useful for bandstructure calculations)

source
DFTK.ExternalFromRealType

External potential from an analytic function V (in Cartesian coordinates). No low-pass filtering is performed.

source
DFTK.FixedBandsType

In 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).

source
DFTK.GPUMethod

Construct a particular GPU architecture by passing the ArrayType

source
DFTK.HartreeType

Hartree 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 = ∑_R 4π δ_R$, integral of G zero on a unit cell).

source
DFTK.KerkerMixingType

Kerker 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)
source
DFTK.KineticType

Kinetic energy:

\[1/2 ∑_n f_n ∫ |∇ψ_n|^2 * {\rm blowup}(-i∇Ψ_n).\]

source
DFTK.KpointType

Discretization 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.

source
DFTK.LazyHcatType

Simple 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).

source
DFTK.LdosModelType

Represents 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.

source
DFTK.LibxcDensitiesMethod
LibxcDensities(
    basis,
    max_derivative::Integer,
    ρ,
    τ
) -> DFTK.LibxcDensities

Compute density in real space and its derivatives starting from ρ

source
DFTK.MagneticType

Magnetic term $A⋅(-i∇)$. It is assumed (but not checked) that $∇⋅A = 0$.

source
DFTK.ModelMethod
Model(system::AbstractSystem; kwargs...)

AtomsBase-compatible Model constructor. Sets structural information (atoms, positions, lattice, n_electrons etc.) from the passed system.

source
DFTK.ModelMethod
Model(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.

source
DFTK.ModelMethod
Model(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.

source
DFTK.NonlocalOperatorType

Nonlocal operator in Fourier space in Kleinman-Bylander format, defined by its projectors P matrix and coupling terms D: Hψ = PDP' ψ.

source
DFTK.NoopOperatorType

Noop operation: don't do anything. Useful for energy terms that don't depend on the orbitals at all (eg nuclei-nuclei interaction).

source
DFTK.PairwisePotentialMethod
PairwisePotential(
    V,
    params;
    max_radius
) -> PairwisePotential

Pairwise 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.

source
DFTK.PlaneWaveBasisType

A 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.

source
DFTK.PlaneWaveBasisMethod

Creates a PlaneWaveBasis using the kinetic energy cutoff Ecut and a k-point grid. By default a MonkhorstPack grid is employed, which can be specified as a MonkhorstPack object or by simply passing a vector of three integers as the kgrid. Optionally kshift allows to specify a shift (0 or 1/2 in each direction). If not specified a grid is generated using kgrid_from_maximal_spacing with a maximal spacing of 2π * 0.022 per Bohr.

source
DFTK.PspHghMethod
PspHgh(path[, identifier, description])

Construct a Hartwigsen, Goedecker, Teter, Hutter separable dual-space Gaussian pseudopotential (1998) from file.

source
DFTK.PspUpfMethod
PspUpf(path[, identifier])

Construct a Unified Pseudopotential Format pseudopotential from file.

Does not support:

  • Fully-realtivistic / spin-orbit pseudos
  • Bare Coulomb / all-electron potentials
  • Semilocal potentials
  • Ultrasoft potentials
  • Projector-augmented wave potentials
  • GIPAW reconstruction data
source
DFTK.RealFourierOperatorType

Linear 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.

source
DFTK.ScfSaveCheckpointsType

Adds checkpointing to a DFTK self-consistent field calculation. The checkpointing file is silently overwritten. Requires the package for writing the output file (usually JLD2) to be loaded.

  • filename: Name of the checkpointing file.
  • compress: Should compression be used on writing (rarely useful)
  • save_ψ: Should the bands also be saved (noteworthy additional cost ... use carefully)
source
DFTK.XcType

Exchange-correlation term, defined by a list of functionals and usually evaluated through libxc.

source
DFTK.χ0MixingType

Generic 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).

source
AbstractFFTs.fft!Method
fft!(
    f_fourier::AbstractArray{T, 3} where T,
    basis::PlaneWaveBasis,
    f_real::AbstractArray{T, 3} where T
) -> Any

In-place version of fft!. NOTE: If kpt is given, not only f_fourier but also f_real is overwritten.

source
AbstractFFTs.fftMethod
fft(
    basis::PlaneWaveBasis{T, VT} where VT<:Real,
    f_real::AbstractArray{U}
) -> Any
fft(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.

source
AbstractFFTs.ifft!Method
ifft!(
    f_real::AbstractArray{T, 3} where T,
    basis::PlaneWaveBasis,
    f_fourier::AbstractArray{T, 3} where T
) -> Any

In-place version of ifft.

source
AbstractFFTs.ifftMethod
ifft(basis::PlaneWaveBasis, f_fourier::AbstractArray) -> Any
ifft(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.

source
AtomsBase.atomic_systemFunction
atomic_system(
    lattice::AbstractMatrix{<:Number},
    atoms::Vector{<:DFTK.Element},
    positions::AbstractVector
) -> AtomsBase.FlexibleSystem
atomic_system(
    lattice::AbstractMatrix{<:Number},
    atoms::Vector{<:DFTK.Element},
    positions::AbstractVector,
    magnetic_moments::AbstractVector
) -> AtomsBase.FlexibleSystem
atomic_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.

source
AtomsBase.periodic_systemFunction
periodic_system(model::Model) -> AtomsBase.FlexibleSystem
periodic_system(
    model::Model,
    magnetic_moments
) -> AtomsBase.FlexibleSystem
periodic_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.

source
Brillouin.KPaths.irrfbz_pathFunction
irrfbz_path(
    model::Model;
    ...
) -> Union{Brillouin.KPaths.KPath{3}, Brillouin.KPaths.KPath{2}, Brillouin.KPaths.KPath{1}}
irrfbz_path(
    model::Model,
    magnetic_moments;
    dim,
    space_group_number
) -> Union{Brillouin.KPaths.KPath{3}, Brillouin.KPaths.KPath{2}, Brillouin.KPaths.KPath{1}}

Extract 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 space_group_number keyword argument.

source
DFTK.CROPFunction
CROP(
    f,
    x0,
    m::Int64,
    max_iter::Int64,
    tol::Real
) -> NamedTuple{(:fixpoint, :converged), <:Tuple{Any, Any}}
CROP(
    f,
    x0,
    m::Int64,
    max_iter::Int64,
    tol::Real,
    warming
) -> NamedTuple{(:fixpoint, :converged), <:Tuple{Any, Any}}

CROP-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.

source
DFTK.G_vectorsMethod
G_vectors(
    basis::PlaneWaveBasis
) -> AbstractArray{StaticArraysCore.SVector{3, Int64}, 3}
G_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.

source
DFTK.G_vectorsMethod
G_vectors(fft_size::Union{Tuple, AbstractVector}) -> Any
G_vectors(fft_size::Tuple)

The wave vectors G in reduced (integer) coordinates for a cubic basis set of given sizes.

source
DFTK.G_vectors_cartMethod
G_vectors_cart(basis::PlaneWaveBasis) -> Any
G_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.

source
DFTK.Gplusk_vectorsMethod
Gplusk_vectors(basis::PlaneWaveBasis, kpt::Kpoint) -> Any
Gplusk_vectors(basis::PlaneWaveBasis, kpt::Kpoint)

The list of $G + k$ vectors, in reduced coordinates.

source
DFTK.Gplusk_vectors_cartMethod
Gplusk_vectors_cart(
    basis::PlaneWaveBasis,
    kpt::Kpoint
) -> Any
Gplusk_vectors_cart(basis::PlaneWaveBasis, kpt::Kpoint)

The list of $G + k$ vectors, in Cartesian coordinates.

source
DFTK.Gplusk_vectors_in_supercellMethod
Gplusk_vectors_in_supercell(
    basis::PlaneWaveBasis,
    basis_supercell::PlaneWaveBasis,
    kpt::Kpoint
) -> Any

Maps all $k+G$ vectors of an given basis as $G$ vectors of the supercell basis, in reduced coordinates.

source
DFTK.HybridMixingMethod
HybridMixing(
;
    εr,
    kTF,
    localization,
    adjust_temperature,
    kwargs...
) -> χ0Mixing

The 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.

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
source
DFTK.IncreaseMixingTemperatureMethod
IncreaseMixingTemperature(
;
    factor,
    above_ρdiff,
    temperature_max
) -> DFTK.var"#callback#688"{DFTK.var"#callback#687#689"{Int64, Float64}}

Increase 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.

source
DFTK.LdosMixingMethod
LdosMixing(; adjust_temperature, kwargs...) -> χ0Mixing

The 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.

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
source
DFTK.ScfAcceptImprovingStepMethod
ScfAcceptImprovingStep(
;
    max_energy_change,
    max_relative_residual
) -> DFTK.var"#accept_step#778"{Float64, Float64}

Accept 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.

source
DFTK.apply_KMethod
apply_K(basis::PlaneWaveBasis, δψ, ψ, ρ, occupation) -> Any
apply_K(basis::PlaneWaveBasis, δψ, ψ, ρ, occupation)

Compute the application of K defined at ψ to δψ. ρ is the density issued from ψ. δψ also generates a δρ, computed with compute_δρ.

source
DFTK.apply_kernelMethod
apply_kernel(
    basis::PlaneWaveBasis,
    δρ;
    RPA,
    kwargs...
) -> Any
apply_kernel(basis::PlaneWaveBasis, δρ; kwargs...)

Computes the potential response to a perturbation δρ in real space, as a 4D (i,j,k,σ) array.

source
DFTK.apply_symopMethod
apply_symop(
    symop::SymOp,
    basis,
    kpoint,
    ψk::AbstractVecOrMat
) -> Tuple{Any, Any}

Apply 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).

source
DFTK.apply_symopMethod
apply_symop(symop::SymOp, basis, ρin; kwargs...) -> Any

Apply a symmetry operation to a density.

source
DFTK.apply_ΩMethod
apply_Ω(δψ, ψ, H::Hamiltonian, Λ) -> Any
apply_Ω(δψ, ψ, 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.

source
DFTK.apply_χ0Method
apply_χ0(
    ham,
    ψ,
    occupation,
    εF,
    eigenvalues,
    δV::AbstractArray{TδV};
    occupation_threshold,
    q,
    kwargs_sternheimer...
) -> Any

Get the density variation δρ corresponding to a potential variation δV.

source
DFTK.attach_pspMethod
attach_psp(
    system::AtomsBase.AbstractSystem,
    pspmap::AbstractDict{Symbol, String}
) -> AtomsBase.FlexibleSystem
attach_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")
source
DFTK.band_data_to_dictMethod
band_data_to_dict(
    band_data::NamedTuple;
    kwargs...
) -> Dict{String, Any}

Convert a band computational result to a dictionary representation. Intended to give a condensed set of results and useful metadata for post processing. See also the todict function for the Model and the PlaneWaveBasis, which are called from this function and the outputs merged. Note, that only the master process returns meaningful data. All other processors still return a dictionary (to simplify code in calling locations), but the data may be dummy.

Some details on the conventions for the returned data:

  • εF: Computed Fermi level (if present in band_data)
  • labels: A mapping of high-symmetry k-Point labels to the index in the kcoords vector of the corresponding k-point.
  • eigenvalues, eigenvalues_error, occupation, residual_norms: (n_bands, n_kpoints, n_spin) arrays of the respective data.
  • n_iter: (n_kpoints, n_spin) array of the number of iterations the diagonalization routine required.
  • kpt_max_n_G: Maximal number of G-vectors used for any k-point.
  • kpt_n_G_vectors: (n_kpoints, n_spin) array, the number of valid G-vectors for each k-point, i.e. the extend along the first axis of ψ where data is valid.
  • kpt_G_vectors: (3, max_n_G, n_kpoints, n_spin) array of the integer (reduced) coordinates of the G-points used for each k-point.
  • ψ: (max_n_G, n_bands, n_kpoints, n_spin) arrays where max_n_G is the maximal number of G-vectors used for any k-point. The data is zero-padded, i.e. for k-points which have less G-vectors than maxnG, then there are tailing zeros.
source
DFTK.build_fft_plans!Method
build_fft_plans!(
    tmp::Array{ComplexF64}
) -> Tuple{FFTW.cFFTWPlan{ComplexF64, -1, true}, FFTW.cFFTWPlan{ComplexF64, -1, false}, Any, Any}

Plan 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.

source
DFTK.build_form_factorsMethod
build_form_factors(
    psp,
    G_plus_k::AbstractArray{StaticArraysCore.SArray{Tuple{3}, TT, 1, 3}, 1}
) -> Matrix{T} where T<:(Complex{_A} where _A)

Build form factors (Fourier transforms of projectors) for an atom centered at 0.

source
DFTK.build_projection_vectorsMethod
build_projection_vectors(
    basis::PlaneWaveBasis{T, VT} where VT<:Real,
    kpt::Kpoint,
    psps,
    psp_positions
) -> Any

Build projection vectors for a atoms array generated by term_nonlocal

\[\begin{aligned} H_{\rm at} &= \sum_{ij} C_{ij} \ket{{\rm proj}_i} \bra{{\rm proj}_j} \\ H_{\rm per} &= \sum_R \sum_{ij} C_{ij} \ket{{\rm proj}_i(x-R)} \bra{{\rm proj}_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} \widehat{\rm proj}_i(k+G') \widehat{\rm proj}_j^*(k+G), \end{aligned}\]

where $\widehat{\rm proj}_i(p) = ∫_{ℝ^3} {\rm proj}_i(r) e^{-ip·r} dr$.

We store $\frac{1}{\sqrt Ω} \widehat{\rm proj}_i(k+G)$ in proj_vectors.

source
DFTK.cell_to_supercellMethod
cell_to_supercell(scfres::NamedTuple) -> NamedTuple

Transpose 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.

source
DFTK.cell_to_supercellMethod
cell_to_supercell(
    basis::PlaneWaveBasis{T, VT} where VT<:Real
) -> PlaneWaveBasis{T, _A, Arch, _B, _C} where {T<:Real, _A<:Real, Arch<:DFTK.AbstractArchitecture, _B<:AbstractArray{StaticArraysCore.SVector{3, Int64}, 3}, _C<:AbstractArray{StaticArraysCore.SVector{3, _A}, 3}}

Construct 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.

source
DFTK.cell_to_supercellMethod
cell_to_supercell(
    ψ,
    basis::PlaneWaveBasis{T<:Real, VT} where VT<:Real,
    basis_supercell::PlaneWaveBasis{T<:Real, VT} where VT<:Real
) -> Any

Re-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.

source
DFTK.cg!Method
cg!(
    x::AbstractArray{T, 1},
    A::LinearMaps.LinearMap{T},
    b::AbstractArray{T, 1};
    precon,
    proj,
    callback,
    tol,
    maxiter,
    miniter
) -> NamedTuple{(:x, :converged, :tol, :residual_norm, :n_iter, :maxiter, :stage), <:Tuple{AbstractVector, Bool, Float64, Any, Int64, Int64, Symbol}}

Implementation of the conjugate gradient method which allows for preconditioning and projection operations along iterations.

source
DFTK.charge_ionicMethod
charge_ionic(el::DFTK.Element) -> Int64

Return the total ionic charge of an atom type (nuclear charge - core electrons)

source
DFTK.compute_amn_kpointMethod
compute_amn_kpoint(
    basis::PlaneWaveBasis,
    kpt,
    ψk,
    projections,
    n_bands
) -> Any

Compute the starting matrix for Wannierization.

Wannierization searches for a unitary matrix $U_{m_n}$. As a starting point for the search, we can provide an initial guess function $g$ for the shape of the Wannier functions, based on what we expect from knowledge of the problem or physical intuition. This starting matrix is called $[A_k]_{m,n}$, and is computed as follows: $[A_k]_{m,n} = \langle ψ_m^k | g^{\text{per}}_n \rangle$. The matrix will be orthonormalized by the chosen Wannier program, we don't need to do so ourselves.

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.

Each projection is a callable object that accepts the basis and some p-points as an argument, and returns the Fourier transform of $g_n$ at the p-points.

source
DFTK.compute_bandsMethod
compute_bands(
    basis_or_scfres,
    kpath::Brillouin.KPaths.KPath;
    kline_density,
    kwargs...
) -> NamedTuple

Compute band data along a specific Brillouin.KPath using a kline_density, the number of $k$-points per inverse bohrs (i.e. overall in units of length).

If not given, the path is determined automatically by inspecting the Model. If you are using spin, you should pass the magnetic_moments as a kwarg to ensure these are taken into account when determining the path.

source
DFTK.compute_bandsMethod
compute_bands(
    scfres::NamedTuple,
    kgrid::DFTK.AbstractKgrid;
    n_bands,
    kwargs...
) -> NamedTuple{(:basis, :ψ, :eigenvalues, :ρ, :εF, :occupation, :diagonalization), <:Tuple{PlaneWaveBasis{T, _A, Arch, _B, _C} where {T<:Real, _A<:Real, Arch<:DFTK.AbstractArchitecture, _B<:AbstractArray{StaticArraysCore.SVector{3, Int64}, 3}, _C<:AbstractArray{StaticArraysCore.SVector{3, _A}, 3}}, Any, Any, Any, Any, Any, Vector{T} where T<:NamedTuple}}

Compute band data starting from SCF results. εF and ρ from the scfres are forwarded to the band computation and n_bands is by default selected as n_bands_scf + 5sqrt(n_bands_scf).

source
DFTK.compute_bandsMethod
compute_bands(
    basis::PlaneWaveBasis,
    kgrid::DFTK.AbstractKgrid;
    n_bands,
    n_extra,
    ρ,
    εF,
    eigensolver,
    tol,
    kwargs...
) -> NamedTuple{(:basis, :ψ, :eigenvalues, :ρ, :εF, :occupation, :diagonalization), <:Tuple{PlaneWaveBasis{T, _A, Arch, _B, _C} where {T<:Real, _A<:Real, Arch<:DFTK.AbstractArchitecture, _B<:AbstractArray{StaticArraysCore.SVector{3, Int64}, 3}, _C<:AbstractArray{StaticArraysCore.SVector{3, _A}, 3}}, Any, Any, Any, Nothing, Nothing, Vector{T} where T<:NamedTuple}}

Compute n_bands eigenvalues and Bloch waves at the k-Points specified by the kgrid. All kwargs not specified below are passed to diagonalize_all_kblocks:

  • kgrid: A custom kgrid to perform the band computation, e.g. a new MonkhorstPack grid.
  • tol The default tolerance for the eigensolver is substantially lower than for SCF computations. Increase if higher accuracy desired.
  • eigensolver: The diagonalisation method to be employed.
source
DFTK.compute_currentMethod
compute_current(
    basis::PlaneWaveBasis,
    ψ,
    occupation
) -> Vector

Computes the probability (not charge) current, $∑_n f_n \Im(ψ_n · ∇ψ_n)$.

source
DFTK.compute_densityMethod
compute_density(
    basis::PlaneWaveBasis{T, VT} where VT<:Real,
    ψ,
    occupation;
    occupation_threshold
) -> AbstractArray{_A, 4} where _A
compute_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.

source
DFTK.compute_dosMethod
compute_dos(
    ε,
    basis,
    eigenvalues;
    smearing,
    temperature
) -> Any

Total density of states at energy ε

source
DFTK.compute_dynmatMethod
compute_dynmat(
    basis::PlaneWaveBasis{T, VT} where VT<:Real,
    ψ,
    occupation;
    q,
    ρ,
    ham,
    εF,
    eigenvalues,
    kwargs...
) -> Any

Compute the dynamical matrix in the form of a $3×n_{\rm atoms}×3×n_{\rm atoms}$ tensor in reduced coordinates.

source
DFTK.compute_fft_sizeMethod
compute_fft_size(
    model::Model{T},
    Ecut;
    ...
) -> Tuple{Vararg{Any, _A}} where _A
compute_fft_size(
    model::Model{T},
    Ecut,
    kgrid;
    ensure_smallprimes,
    algorithm,
    factors,
    kwargs...
) -> Tuple{Vararg{Any, _A}} where _A

Determine 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

source
DFTK.compute_forcesMethod
compute_forces(scfres) -> Any

Compute 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.

source
DFTK.compute_forces_cartMethod
compute_forces_cart(
    basis::PlaneWaveBasis,
    ψ,
    occupation;
    kwargs...
) -> Any

Compute 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.

source
DFTK.compute_inverse_latticeMethod
compute_inverse_lattice(lattice::AbstractArray{T, 2}) -> Any

Compute the inverse of the lattice. Takes special care of 1D or 2D cases.

source
DFTK.compute_kernelMethod
compute_kernel(
    basis::PlaneWaveBasis{T, VT} where VT<:Real;
    kwargs...
) -> Any
compute_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)\]

source
DFTK.compute_ldosMethod
compute_ldos(
    ε,
    basis::PlaneWaveBasis{T, VT} where VT<:Real,
    eigenvalues,
    ψ;
    smearing,
    temperature,
    weight_threshold
) -> AbstractArray{_A, 4} where _A

Local density of states, in real space. weight_threshold is a threshold to screen away small contributions to the LDOS.

source
DFTK.compute_occupationMethod
compute_occupation(
    basis::PlaneWaveBasis{T, VT} where VT<:Real,
    eigenvalues::AbstractVector,
    εF::Number;
    temperature,
    smearing
) -> NamedTuple{(:occupation, :εF), <:Tuple{Any, Number}}

Compute occupation given eigenvalues and Fermi level

source
DFTK.compute_occupationMethod
compute_occupation(
    basis::PlaneWaveBasis{T, VT} where VT<:Real,
    eigenvalues::AbstractVector;
    ...
) -> NamedTuple{(:occupation, :εF), <:Tuple{Any, Number}}
compute_occupation(
    basis::PlaneWaveBasis{T, VT} where VT<:Real,
    eigenvalues::AbstractVector,
    fermialg::AbstractFermiAlgorithm;
    tol_n_elec,
    temperature,
    smearing
) -> NamedTuple{(:occupation, :εF), <:Tuple{Any, Number}}

Compute 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.

source
DFTK.compute_recip_latticeMethod
compute_recip_lattice(lattice::AbstractArray{T, 2}) -> Any

Compute 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.

source
DFTK.compute_stresses_cartMethod
compute_stresses_cart(scfres) -> Any

Compute the stresses of an obtained SCF solution. The stress tensor is given by

\[\left( \begin{array}{ccc} σ_{xx} σ_{xy} σ_{xz} \\ σ_{yx} σ_{yy} σ_{yz} \\ σ_{zx} σ_{zy} σ_{zz} \end{array} \right) = \frac{1}{|Ω|} \left. \frac{dE[ (I+ϵ) * L]}{dM}\right|_{ϵ=0}\]

where $ϵ$ is the strain. See O. Nielsen, R. Martin Phys. Rev. B. 32, 3792 (1985) for details. In Voigt notation one would use the vector $[σ_{xx} σ_{yy} σ_{zz} σ_{zy} σ_{zx} σ_{yx}]$.

source
DFTK.compute_transfer_matrixMethod
compute_transfer_matrix(
    basis_in::PlaneWaveBasis,
    kpt_in::Kpoint,
    basis_out::PlaneWaveBasis,
    kpt_out::Kpoint
) -> Any

Return a sparse matrix that maps quantities given on basis_in and kpt_in to quantities on basis_out and kpt_out.

source
DFTK.compute_transfer_matrixMethod
compute_transfer_matrix(
    basis_in::PlaneWaveBasis,
    basis_out::PlaneWaveBasis
) -> Vector

Return 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.

source
DFTK.compute_unit_cell_volumeMethod
compute_unit_cell_volume(lattice) -> Any

Compute unit cell volume volume. In case of 1D or 2D case, the volume is the length/surface.

source
DFTK.compute_δHψ_αsMethod
compute_δHψ_αs(basis::PlaneWaveBasis, ψ, α, s, q) -> Any

Assemble the right-hand side term for the Sternheimer equation for all relevant quantities: Compute the perturbation of the Hamiltonian with respect to a variation of the local potential produced by a displacement of the atom s in the direction α.

source
DFTK.compute_δocc!Method
compute_δocc!(
    δocc,
    basis::PlaneWaveBasis{T, VT} where VT<:Real,
    ψ,
    εF,
    ε,
    δHψ
) -> NamedTuple{(:δocc, :δεF), <:Tuple{Any, Any}}

Compute 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.

source
DFTK.compute_δψ!Method
compute_δψ!(
    δψ,
    basis::PlaneWaveBasis{T, VT} where VT<:Real,
    H,
    ψ,
    εF,
    ε,
    δHψ;
    ...
)
compute_δψ!(
    δψ,
    basis::PlaneWaveBasis{T, VT} where VT<:Real,
    H,
    ψ,
    εF,
    ε,
    δHψ,
    ε_minus_q;
    ψ_extra,
    q,
    kwargs_sternheimer...
)

Perform 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.

source
DFTK.compute_χ0Method
compute_χ0(ham; temperature) -> Any

Compute 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)\]

source
DFTK.count_n_projMethod
count_n_proj(psps, psp_positions) -> Any
count_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.

source
DFTK.count_n_projMethod
count_n_proj(psp::DFTK.NormConservingPsp, l::Integer) -> Any
count_n_proj(psp, l)

Number of projector functions for angular momentum l, including angular parts from -m:m.

source
DFTK.count_n_projMethod
count_n_proj(psp::DFTK.NormConservingPsp) -> Int64
count_n_proj(psp)

Number of projector functions for all angular momenta up to psp.lmax, including angular parts from -m:m.

source
DFTK.count_n_proj_radialMethod
count_n_proj_radial(
    psp::DFTK.NormConservingPsp,
    l::Integer
) -> Any
count_n_proj_radial(psp, l)

Number of projector radial functions at angular momentum l.

source
DFTK.count_n_proj_radialMethod
count_n_proj_radial(psp::DFTK.NormConservingPsp) -> Int64
count_n_proj_radial(psp)

Number of projector radial functions at all angular momenta up to psp.lmax.

source
DFTK.create_supercellMethod
create_supercell(
    lattice,
    atoms,
    positions,
    supercell_size
) -> NamedTuple{(:lattice, :atoms, :positions), <:Tuple{Any, Any, Any}}

Construct a supercell of size supercell_size from a unit cell described by its lattice, atoms and their positions.

source
DFTK.default_fermialgMethod
default_fermialg(
    _::DFTK.Smearing.SmearingFunction
) -> FermiBisection

Default selection of a Fermi level determination algorithm

source
DFTK.default_symmetriesMethod
default_symmetries(
    lattice,
    atoms,
    positions,
    magnetic_moments,
    spin_polarization,
    terms;
    tol_symmetry
) -> Union{Vector{SymOp{Bool}}, Vector{SymOp{Float64}}}

Default logic to determine the symmetry operations to be used in the model.

source
DFTK.default_wannier_centersMethod
default_wannier_centers(n_wannier) -> Any

Default random Gaussian guess for maximally-localised wannier functions generated in reduced coordinates.

source
DFTK.diagonalize_all_kblocksMethod
diagonalize_all_kblocks(
    eigensolver,
    ham::Hamiltonian,
    nev_per_kpoint::Int64;
    ψguess,
    prec_type,
    interpolate_kpoints,
    tol,
    miniter,
    maxiter,
    n_conv_check
) -> NamedTuple{(:λ, :X, :residual_norms, :n_iter, :converged, :n_matvec), <:Tuple{Vector, Vector, Vector, Vector, Union{Missing, Bool}, Any}}

Function 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)

source
DFTK.diameterMethod
diameter(lattice::AbstractMatrix) -> Any

Compute the diameter of the unit cell

source
DFTK.direct_minimizationMethod
direct_minimization(
    basis::PlaneWaveBasis{T, VT} where VT<:Real;
    ψ,
    tol,
    is_converged,
    maxiter,
    prec_type,
    callback,
    optim_method,
    linesearch,
    kwargs...
) -> Any

Computes the ground state by direct minimization. kwargs... are passed to Optim.Options() and optim_method selects the optim approach which is employed.

source
DFTK.disable_threadingMethod
disable_threading() -> Union{Nothing, Bool}

Convenience function to disable all threading in DFTK and assert that Julia threading is off as well.

source
DFTK.divergence_realMethod
divergence_real(operand, basis) -> Any

Compute 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.

source
DFTK.energy_forces_ewaldMethod
energy_forces_ewald(
    lattice::AbstractArray{T},
    charges::AbstractArray,
    positions;
    kwargs...
) -> NamedTuple{(:energy, :forces), <:Tuple{Any, Any}}

Standard computation of energy and forces.

source
DFTK.energy_forces_ewaldMethod
energy_forces_ewald(
    lattice::AbstractArray{T},
    charges,
    positions,
    q,
    ph_disp;
    kwargs...
) -> NamedTuple{(:energy, :forces), <:Tuple{Any, Any}}

Computation for phonons; required to build the dynamical matrix.

source
DFTK.energy_forces_ewaldMethod
energy_forces_ewald(
    S,
    lattice::AbstractArray{T},
    charges,
    positions,
    q,
    ph_disp;
    η
) -> NamedTuple{(:energy, :forces), <:Tuple{Any, Any}}

Compute 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.

source
DFTK.energy_forces_pairwiseMethod
energy_forces_pairwise(
    lattice::AbstractArray{T},
    symbols,
    positions,
    V,
    params;
    kwargs...
) -> NamedTuple{(:energy, :forces), <:Tuple{Any, Any}}

Standard computation of energy and forces.

source
DFTK.energy_forces_pairwiseMethod
energy_forces_pairwise(
    lattice::AbstractArray{T},
    symbols,
    positions,
    V,
    params,
    q,
    ph_disp;
    kwargs...
) -> NamedTuple{(:energy, :forces), <:Tuple{Any, Any}}

Computation for phonons; required to build the dynamical matrix.

source
DFTK.energy_forces_pairwiseMethod
energy_forces_pairwise(
    S,
    lattice::AbstractArray{T},
    symbols,
    positions,
    V,
    params,
    q,
    ph_disp;
    max_radius
) -> NamedTuple{(:energy, :forces), <:Tuple{Any, Any}}

Compute 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.

source
DFTK.energy_psp_correctionMethod
energy_psp_correction(
    lattice::AbstractArray{T, 2},
    atoms,
    atom_groups
) -> Any

Compute the correction term for properly modelling the interaction of the pseudopotential core with the compensating background charge induced by the Ewald term.

source
DFTK.enforce_real!Method
enforce_real!(
    fourier_coeffs,
    basis::PlaneWaveBasis
) -> AbstractArray

Ensure 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.

source
DFTK.estimate_integer_lattice_boundsMethod
estimate_integer_lattice_bounds(
    M::AbstractArray{T, 2},
    δ;
    ...
) -> Vector
estimate_integer_lattice_bounds(
    M::AbstractArray{T, 2},
    δ,
    shift;
    tol
) -> Vector

Estimate 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.

source
DFTK.eval_psp_density_core_fourierMethod
eval_psp_density_core_fourier(
    _::DFTK.NormConservingPsp,
    _::Real
) -> Any
eval_psp_density_core_fourier(psp, p)

Evaluate the atomic core charge density in reciprocal space:

\[\begin{aligned} ρ_{\rm core}(p) &= ∫_{ℝ^3} ρ_{\rm core}(r) e^{-ip·r} dr \\ &= 4π ∫_{ℝ_+} ρ_{\rm core}(r) \frac{\sin(p·r)}{ρ·r} r^2 dr. \end{aligned}\]

source
DFTK.eval_psp_density_core_realMethod
eval_psp_density_core_real(
    _::DFTK.NormConservingPsp,
    _::Real
) -> Any
eval_psp_density_core_real(psp, r)

Evaluate the atomic core charge density in real space.

source
DFTK.eval_psp_density_valence_fourierMethod
eval_psp_density_valence_fourier(
    psp::DFTK.NormConservingPsp,
    p::AbstractVector
) -> Any
eval_psp_density_valence_fourier(psp, p)

Evaluate the atomic valence charge density in reciprocal space:

\[\begin{aligned} ρ_{\rm val}(p) &= ∫_{ℝ^3} ρ_{\rm val}(r) e^{-ip·r} dr \\ &= 4π ∫_{ℝ_+} ρ_{\rm val}(r) \frac{\sin(p·r)}{ρ·r} r^2 dr. \end{aligned}\]

source
DFTK.eval_psp_density_valence_realMethod
eval_psp_density_valence_real(
    psp::DFTK.NormConservingPsp,
    r::AbstractVector
) -> Any
eval_psp_density_valence_real(psp, r)

Evaluate the atomic valence charge density in real space.

source
DFTK.eval_psp_energy_correctionFunction
eval_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 $p \to 0$, using the same correction as in the Fourier-transform of the local potential:

\[\lim_{p \to 0} 4π N_{\rm elec} ∫_{ℝ_+} (V(r) - C(r)) \frac{\sin(p·r)}{p·r} r^2 dr + F[C(r)] = 4π N_{\rm elec} ∫_{ℝ_+} (V(r) + Z/r) r^2 dr.\]

source
DFTK.eval_psp_local_fourierMethod
eval_psp_local_fourier(
    psp::DFTK.NormConservingPsp,
    p::AbstractVector
) -> Any
eval_psp_local_fourier(psp, p)

Evaluate the local part of the pseudopotential in reciprocal space:

\[\begin{aligned} V_{\rm loc}(p) &= ∫_{ℝ^3} V_{\rm loc}(r) e^{-ip·r} dr \\ &= 4π ∫_{ℝ_+} V_{\rm loc}(r) \frac{\sin(p·r)}{p} 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}(p) &= ∫_{ℝ^3} (V_{\rm loc}(r) - C(r)) e^{-ip·r} dr + F[C(r)] \\ &= 4π ∫_{ℝ_+} (V_{\rm loc}(r) + Z/r) \frac{\sin(p·r)}{p·r} r^2 dr - Z/p^2. \end{aligned}\]

source
DFTK.eval_psp_local_realMethod
eval_psp_local_real(
    psp::DFTK.NormConservingPsp,
    r::AbstractVector
) -> Any
eval_psp_local_real(psp, r)

Evaluate the local part of the pseudopotential in real space.

source
DFTK.eval_psp_projector_fourierMethod
eval_psp_projector_fourier(
    psp::DFTK.NormConservingPsp,
    p::AbstractVector
)
eval_psp_projector_fourier(psp, i, l, p)

Evaluate the radial part of the i-th projector for angular momentum l at the reciprocal vector with modulus p:

\[\begin{aligned} {\rm proj}(p) &= ∫_{ℝ^3} {\rm proj}_{il}(r) e^{-ip·r} dr \\ &= 4π ∫_{ℝ_+} r^2 {\rm proj}_{il}(r) j_l(p·r) dr. \end{aligned}\]

source
DFTK.eval_psp_projector_realMethod
eval_psp_projector_real(
    psp::DFTK.NormConservingPsp,
    i,
    l,
    r::AbstractVector
) -> Any
eval_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.

source
DFTK.filled_occupationMethod
filled_occupation(model) -> Int64

Maximal occupation of a state (2 for non-spin-polarized electrons, 1 otherwise).

source
DFTK.find_equivalent_kptMethod
find_equivalent_kpt(
    basis::PlaneWaveBasis{T, VT} where VT<:Real,
    kcoord,
    spin;
    tol
) -> NamedTuple{(:index, :ΔG), <:Tuple{Int64, Any}}

Find the equivalent index of the coordinate kcoord ∈ ℝ³ in a list kcoords ∈ [-½, ½)³. ΔG is the vector of ℤ³ such that kcoords[index] = kcoord + ΔG.

source
DFTK.gather_kptsMethod
gather_kpts(
    basis::PlaneWaveBasis
) -> Union{Nothing, PlaneWaveBasis}

Gather 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 for debugging or to extract data for serialisation to disk.

source
DFTK.gather_kpts_block!Method
gather_kpts_block!(
    dest,
    basis::PlaneWaveBasis,
    kdata::AbstractArray{A, 1}
) -> Any

Gather the distributed data of a quantity depending on k-Points on the master process and save it in dest as a dense (size(kdata[1])..., n_kpoints) array. On the other (non-master) processes nothing is returned.

source
DFTK.guess_densityFunction
guess_density(
    basis::PlaneWaveBasis,
    method::AtomicDensity;
    ...
) -> Any
guess_density(
    basis::PlaneWaveBasis,
    method::AtomicDensity,
    magnetic_moments;
    n_electrons
) -> Any
guess_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$.

source
DFTK.hankelMethod
hankel(
    r::AbstractVector,
    r2_f::AbstractVector,
    l::Integer,
    p::Real
) -> Any
hankel(r, r2_f, l, p)

Compute the Hankel transform

\[ H[f] = 4\pi \int_0^\infty r f(r) j_l(p·r) r dr.\]

The integration is performed by trapezoidal quadrature, and the function takes as input the radial grid r, the precomputed quantity r²f(r) r2_f, angular momentum / spherical bessel order l, and the Hankel coordinate p.

source
DFTK.has_core_densityMethod
has_core_density(_::DFTK.Element) -> Any

Check presence of model core charge density (non-linear core correction).

source
DFTK.index_G_vectorsMethod
index_G_vectors(
    fft_size::Tuple,
    G::AbstractVector{<:Integer}
) -> Any

Return 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.

source
DFTK.interpolate_densityMethod
interpolate_density(
    ρ_in::AbstractArray{T, 4},
    basis_in::PlaneWaveBasis,
    basis_out::PlaneWaveBasis
) -> AbstractArray{T, 4} where T

Interpolate a density 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.

source
DFTK.interpolate_densityMethod
interpolate_density(
    ρ_in::AbstractArray{T, 4},
    grid_in::Tuple{T, T, T} where T,
    grid_out::Tuple{T, T, T} where T,
    lattice_in,
    lattice_out
) -> AbstractArray{T, 4} where T

Interpolate a density in real space from one FFT grid to another, where lattice_in and lattice_out may be supercells of each other.

source
DFTK.interpolate_densityMethod
interpolate_density(
    ρ_in::AbstractArray{T, 4},
    grid_out::Tuple{T, T, T} where T
) -> AbstractArray{T, 4} where T

Interpolate a density in real space from one FFT grid to another. Assumes the lattice is unchanged.

source
DFTK.interpolate_kpointMethod
interpolate_kpoint(
    data_in::AbstractVecOrMat,
    basis_in::PlaneWaveBasis,
    kpoint_in::Kpoint,
    basis_out::PlaneWaveBasis,
    kpoint_out::Kpoint
) -> Any

Interpolate some data from one $k$-point to another. The interpolation is fast, but not necessarily exact. Intended only to construct guesses for iterative solvers.

source
DFTK.irfftMethod
irfft(
    basis::PlaneWaveBasis{T, VT} where VT<:Real,
    f_fourier::AbstractArray
) -> Any

Perform a real valued iFFT; see ifft. Note that this function silently drops the imaginary part.

source
DFTK.irreducible_kcoordsMethod
irreducible_kcoords(
    kgrid::MonkhorstPack,
    symmetries::AbstractVector{<:SymOp};
    check_symmetry
) -> Union{@NamedTuple{kcoords::Vector{StaticArraysCore.SVector{3, Rational{Int64}}}, kweights::Vector{Float64}}, @NamedTuple{kcoords::Vector{StaticArraysCore.SVector{3, Float64}}, kweights::Vector{Float64}}}

Construct the irreducible wedge given the crystal symmetries. Returns the list of k-point coordinates and the associated weights.

source
DFTK.is_metalMethod
is_metal(eigenvalues, εF; tol) -> Bool
is_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.

source
DFTK.k_to_equivalent_kpq_permutationMethod
k_to_equivalent_kpq_permutation(
    basis::PlaneWaveBasis,
    qcoord
) -> Vector

Return the indices of the kpoints shifted by q. That is for each kpoint of the basis: kpoints[ik].coordinate + q is equivalent to kpoints[indices[ik]].coordinate.

source
DFTK.kgrid_from_maximal_spacingMethod
kgrid_from_maximal_spacing(
    lattice,
    spacing;
    kshift
) -> Union{MonkhorstPack, Vector{Int64}}

Build a MonkhorstPack grid to ensure kpoints are at most this spacing apart (in inverse Bohrs). A reasonable spacing is 0.13 inverse Bohrs (around $2π * 0.04 \, \text{Å}^{-1}$).

source
DFTK.kgrid_from_minimal_n_kpointsMethod
kgrid_from_minimal_n_kpoints(
    lattice,
    n_kpoints::Integer;
    kshift
) -> Union{MonkhorstPack, Vector{Int64}}

Selects a MonkhorstPack grid 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.

source
DFTK.kpath_get_kcoordsMethod
kpath_get_kcoords(
    kinter::Brillouin.KPaths.KPathInterpolant{D}
) -> Vector

Return kpoint coordinates in reduced coordinates

source
DFTK.kpq_equivalent_blochwave_to_kpqMethod
kpq_equivalent_blochwave_to_kpq(
    basis,
    kpt,
    q,
    ψk_plus_q_equivalent
) -> NamedTuple{(:kpt, :ψk), <:Tuple{Kpoint, Any}}

Create the Fourier expansion of $ψ_{k+q}$ from $ψ_{[k+q]}$, where $[k+q]$ is in basis.kpoints. while $k+q$ may or may not be inside.

If $ΔG ≔ [k+q] - (k+q)$, then we have that

\[ ∑_G \hat{u}_{[k+q]}(G) e^{i(k+q+G)·r} = ∑_{G'} \hat{u}_{k+q}(G'-ΔG) e^{i(k+q+ΔG+G')·r},\]

hence

\[ u_{k+q}(G) = u_{[k+q]}(G + ΔG).\]

source
DFTK.krange_spinMethod
krange_spin(basis::PlaneWaveBasis, spin::Integer) -> Any

Return the index range of $k$-points that have a particular spin component.

source
DFTK.kwargs_scf_checkpointsMethod
kwargs_scf_checkpoints(
    basis::DFTK.AbstractBasis;
    filename,
    callback,
    diagtolalg,
    ρ,
    ψ,
    save_ψ,
    kwargs...
) -> NamedTuple{(:callback, :diagtolalg, :ψ, :ρ), <:Tuple{ComposedFunction{ScfDefaultCallback, ScfSaveCheckpoints}, AdaptiveDiagtol, Any, Any}}

Transparently handle checkpointing by either returning kwargs for self_consistent_field, which start checkpointing (if no checkpoint file is present) or that continue a checkpointed run (if a checkpoint file can be loaded). filename is the location where the checkpoint is saved, save_ψ determines whether orbitals are saved in the checkpoint as well. The latter is discouraged, since generally slow.

source
DFTK.list_pspFunction
list_psp(; ...) -> Any
list_psp(element; family, functional, core) -> Any
list_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.

source
DFTK.load_pspMethod
load_psp(
    key::AbstractString;
    kwargs...
) -> Union{PspHgh{Float64}, PspUpf{_A, Interpolations.Extrapolation{T, 1, ITPT, IT, Interpolations.Throw{Nothing}}} where {_A, T, ITPT, IT}}

Load 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.

source
DFTK.load_scfresFunction
load_scfres(filename::AbstractString; ...) -> Any
load_scfres(
    filename::AbstractString,
    basis;
    skip_hamiltonian,
    strict
) -> Any

Load back an scfres, which has previously been stored with save_scfres. Note the warning in save_scfres.

If basis is nothing, the basis is also loaded and reconstructed from the file, in which case architecture=CPU(). If a basis is passed, this one is used, which can be used to continue computation on a slightly different model or to avoid the cost of rebuilding the basis. If the stored basis and the passed basis are inconsistent (e.g. different FFT size, Ecut, k-points etc.) the load_scfres will error out.

By default the energies and ham (Hamiltonian object) are recomputed. To avoid this, set skip_hamiltonian=true. On errors the routine exits unless strict=false in which case it tries to recover from the file as much data as it can, but then the resulting scfres might not be fully consistent.

No compatibility guarantees

No guarantees are made with respect to this function at this point. It may change incompatibly between DFTK versions (including patch versions) or stop working / be removed in the future.

source
DFTK.model_DFTMethod
model_DFT(
    lattice::AbstractMatrix,
    atoms::Vector{<:DFTK.Element},
    positions::Vector{<:AbstractVector},
    xc::Xc;
    extra_terms,
    kwargs...
) -> Model

Build a DFT model from the specified atoms, with the specified functionals.

source
DFTK.model_atomicMethod
model_atomic(
    lattice::AbstractMatrix,
    atoms::Vector{<:DFTK.Element},
    positions::Vector{<:AbstractVector};
    extra_terms,
    kinetic_blowup,
    kwargs...
) -> Model

Convenience constructor, which builds a standard atomic (kinetic + atomic potential) model. Use extra_terms to add additional terms.

source
DFTK.mpi_nprocsFunction
mpi_nprocs() -> Int64
mpi_nprocs(comm) -> Int64

Number of processors used in MPI. Can be called without ensuring initialization.

source
DFTK.multiply_ψ_by_blochwaveMethod
multiply_ψ_by_blochwave(
    basis::PlaneWaveBasis,
    ψ,
    f_real,
    q
) -> Any
multiply_ψ_by_blochwave(basis::PlaneWaveBasis, ψ, f_real, q)

Return the Fourier coefficients for the Bloch waves $f^{\rm real}_{q} ψ_{k-q}$ in an element of basis.kpoints equivalent to $k-q$.

source
DFTK.newtonMethod
newton(
    basis::PlaneWaveBasis{T, VT} where VT<:Real,
    ψ0;
    tol,
    tol_cg,
    maxiter,
    callback,
    is_converged
) -> NamedTuple{(:ham, :basis, :energies, :converged, :ρ, :eigenvalues, :occupation, :εF, :n_iter, :ψ, :stage, :algorithm, :runtime_ns), <:Tuple{Hamiltonian, PlaneWaveBasis, Energies, Any, AbstractArray{_A, 4} where _A, Vector{Any}, Vector, Nothing, Int64, Any, Symbol, String, UInt64}}
newton(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.

source
DFTK.next_compatible_fft_sizeMethod
next_compatible_fft_size(
    size::Int64;
    smallprimes,
    factors
) -> Int64

Find 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.

source
DFTK.next_densityFunction
next_density(
    ham::Hamiltonian;
    ...
) -> NamedTuple{(:ψ, :eigenvalues, :occupation, :εF, :ρout, :diagonalization, :n_bands_converge, :occupation_threshold), <:Tuple{Vector, Vector, Vector, Number, AbstractArray{_A, 4} where _A, NamedTuple{(:λ, :X, :residual_norms, :n_iter, :converged, :n_matvec), <:Tuple{Vector, Vector, Vector, Vector, Union{Missing, Bool}, Any}}, Int64, Float64}}
next_density(
    ham::Hamiltonian,
    nbandsalg::DFTK.NbandsAlgorithm;
    ...
) -> NamedTuple{(:ψ, :eigenvalues, :occupation, :εF, :ρout, :diagonalization, :n_bands_converge, :occupation_threshold), <:Tuple{Vector, Vector, Vector, Number, AbstractArray{_A, 4} where _A, NamedTuple{(:λ, :X, :residual_norms, :n_iter, :converged, :n_matvec), <:Tuple{Vector, Vector, Vector, Vector, Union{Missing, Bool}, Any}}, Int64, Any}}
next_density(
    ham::Hamiltonian,
    nbandsalg::DFTK.NbandsAlgorithm,
    fermialg::AbstractFermiAlgorithm;
    eigensolver,
    ψ,
    eigenvalues,
    occupation,
    kwargs...
) -> NamedTuple{(:ψ, :eigenvalues, :occupation, :εF, :ρout, :diagonalization, :n_bands_converge, :occupation_threshold), <:Tuple{Vector, Vector, Vector, Number, AbstractArray{_A, 4} where _A, NamedTuple{(:λ, :X, :residual_norms, :n_iter, :converged, :n_matvec), <:Tuple{Vector, Vector, Vector, Vector, Union{Missing, Bool}, Any}}, Int64, Any}}

Obtain 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.

source
DFTK.norm_cplxMethod
norm_cplx(x) -> Any

Complex-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.

source
DFTK.overlap_Mmn_k_kpbMethod
overlap_Mmn_k_kpb(
    basis::PlaneWaveBasis,
    ψ,
    ik,
    ik_plus_b,
    G_shift,
    n_bands
) -> Any

Computes the matrix $[M^{k,b}]_{m,n} = \langle u_{m,k} | u_{n,k+b} \rangle$ for given k.

G_shift is the "shifting" vector, correction due to the periodicity conditions imposed on $k \to ψ_k$. It is non zero if k_plus_b 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}$.

source
DFTK.pcut_psp_localMethod
pcut_psp_local(psp::PspHgh{T}) -> Any

Estimate an upper bound for the argument p after which abs(eval_psp_local_fourier(psp, p)) is a strictly decreasing function.

source
DFTK.pcut_psp_projectorMethod
pcut_psp_projector(psp::PspHgh{T}, i, l) -> Any

Estimate an upper bound for the argument p after which eval_psp_projector_fourier(psp, p) is a strictly decreasing function.

source
DFTK.phonon_modesMethod
phonon_modes(
    basis::PlaneWaveBasis{T, VT} where VT<:Real,
    dynmat
) -> NamedTuple{(:frequencies, :vectors), <:Tuple{Any, Any}}

Solve the eigenproblem for a dynamical matrix: returns the frequencies and eigenvectors (vectors).

source
DFTK.plot_bandstructureFunction

Compute and plot the band structure. Kwargs are like in compute_bands. Requires Plots.jl to be loaded to be defined and working properly. 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).

source
DFTK.plot_dosFunction

Plot the density of states over a reasonable range. Requires to load Plots.jl beforehand.

source
DFTK.psp_local_polynomialFunction
psp_local_polynomial(T, psp::PspHgh) -> Any
psp_local_polynomial(T, psp::PspHgh, t) -> Any

The 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} p$ and Q is a polynomial of at most degree 8. This function returns Q.

source
DFTK.psp_projector_polynomialFunction
psp_projector_polynomial(T, psp::PspHgh, i, l) -> Any
psp_projector_polynomial(T, psp::PspHgh, i, l, t) -> Any

The 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 p$ and Q is a polynomial. This function returns Q.

source
DFTK.r_vectorsMethod
r_vectors(
    basis::PlaneWaveBasis
) -> AbstractArray{StaticArraysCore.SVector{3, VT}, 3} where VT<:Real
r_vectors(basis::PlaneWaveBasis)

The list of $r$ vectors, in reduced coordinates. By convention, this is in [0,1)^3.

source
DFTK.r_vectors_cartMethod
r_vectors_cart(basis::PlaneWaveBasis) -> Any
r_vectors_cart(basis::PlaneWaveBasis)

The list of $r$ vectors, in Cartesian coordinates.

source
DFTK.radial_hydrogenicMethod
radial_hydrogenic(
    r::AbstractArray{T<:Real, 1},
    n::Integer
) -> Any
radial_hydrogenic(
    r::AbstractArray{T<:Real, 1},
    n::Integer,
    α::Real
) -> Any

Radial functions from solutions of Hydrogenic Schrödinger equation. Same as Wannier90 user guide Table 3.3.

Arguments

  • r: radial grid
  • n: principal quantum number
  • α: diffusivity, $\frac{Z}{a}$ where $Z$ is the atomic number and $a$ is the Bohr radius.
source
DFTK.random_densityMethod
random_density(
    basis::PlaneWaveBasis{T, VT} where VT<:Real,
    n_electrons::Integer
) -> Any

Build a random charge density normalized to the provided number of electrons.

source
DFTK.read_w90_nnkpMethod
read_w90_nnkp(
    fileprefix::String
) -> @NamedTuple{nntot::Int64, nnkpts::Vector{@NamedTuple{ik::Int64, ik_plus_b::Int64, G_shift::Vector{Int64}}}}

Read the .nnkp file provided by the preprocessing routine of Wannier90 (i.e. "wannier90.x -pp prefix") Returns:

  1. the array 'nnkpts' of k points, their respective nearest neighbors and associated shifing vectors (non zero if the neighbor is located in another cell).
  2. the number 'nntot' of neighbors per k point.

TODO: add the possibility to exclude bands

source
DFTK.reducible_kcoordsMethod
reducible_kcoords(
    kgrid::MonkhorstPack
) -> @NamedTuple{kcoords::Vector{StaticArraysCore.SVector{3, Rational{Int64}}}}

Construct the coordinates of the k-points in a (shifted) Monkhorst-Pack grid

source
DFTK.run_wannier90Function

Wannerize 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. Random Gaussians are used as guesses by default, can be changed using the projections kwarg. All keyword arguments supported by Wannier90 for the disentanglement may be added as keyword arguments. The function returns the fileprefix.

Experimental feature

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.

source
DFTK.save_bandsMethod
save_bands(
    filename::AbstractString,
    band_data::NamedTuple;
    save_ψ
)

Write the computed bands to a file. On all processes, but the master one the filename is ignored. save_ψ determines whether the wavefunction is also saved or not. Note that this function can be both used on the results of compute_bands and self_consistent_field.

Changes to data format reserved

No guarantees are made with respect to the format of the keys at this point. We may change this incompatibly between DFTK versions (including patch versions). In particular changes with respect to the ψ structure are planned.

source
DFTK.save_scfresMethod
save_scfres(
    filename::AbstractString,
    scfres::NamedTuple;
    save_ψ,
    extra_data,
    compress,
    save_ρ
) -> Any

Save an scfres obtained from self_consistent_field to a file. On all processes but the master one the filename is ignored. 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. Note that this file is also a valid HDF5 file, which can thus similarly be read by external non-Julia libraries such as h5py or similar. 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, optionally bands and some metadata.
  • json: A JSON file with basic information about the SCF run. Stores for example the number of iterations, occupations, some information about the basis, eigenvalues, Fermi level etc.

Keyword arguments:

  • save_ψ: Save the orbitals as well (may lead to larger files). This is the default for jld2, but false for all other formats, where this is considerably more expensive.
  • save_ρ: Save the density as well (may lead to larger files). This is the default for all but json.
  • extra_data: Additional data to place into the file. The data is just copied like fp["key"] = value, where fp is a JLD2.JLDFile, WriteVTK.vtk_grid and so on.
  • compress: Apply compression to array data. Requires the CodecZlib package to be available.
Changes to data format reserved

No guarantees are made with respect to the format of the keys at this point. We may change this incompatibly between DFTK versions (including patch versions). In particular changes with respect to the ψ structure are planned.

source
DFTK.scatter_kpts_blockMethod
scatter_kpts_block(
    basis::PlaneWaveBasis,
    data::Union{Nothing, AbstractArray}
) -> Any

Scatter the data of a quantity depending on k-Points from the master process to the child processes and return it as a Vector{Array}, where the outer vector is a list over all k-points. On non-master processes nothing may be passed.

source
DFTK.scf_anderson_solverFunction
scf_anderson_solver(
;
    ...
) -> DFTK.var"#anderson#695"{DFTK.var"#anderson#694#696"{Base.Pairs{Symbol, Union{}, Tuple{}, @NamedTuple{}}, Int64}}
scf_anderson_solver(
    m;
    kwargs...
) -> DFTK.var"#anderson#695"{DFTK.var"#anderson#694#696"{Base.Pairs{Symbol, Union{}, Tuple{}, @NamedTuple{}}, _A}} where _A

Create a simple anderson-accelerated SCF solver. m specifies the number of steps to keep the history of.

source
DFTK.scf_damping_quadratic_modelMethod
scf_damping_quadratic_model(
    info,
    info_next;
    modeltol
) -> NamedTuple{(:α, :relerror), <:Tuple{Any, Any}}

Use 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.

source
DFTK.scf_damping_solverFunction
scf_damping_solver(

) -> DFTK.var"#fp_solver#691"{DFTK.var"#fp_solver#690#692"}
scf_damping_solver(
    β
) -> DFTK.var"#fp_solver#691"{DFTK.var"#fp_solver#690#692"}

Create a damped SCF solver updating the density as x = β * x_new + (1 - β) * x

source
DFTK.scfres_to_dictMethod
scfres_to_dict(
    scfres::NamedTuple;
    kwargs...
) -> Dict{String, Any}

Convert an scfres to a dictionary representation. Intended to give a condensed set of results and useful metadata for post processing. See also the todict function for the Model and the PlaneWaveBasis as well as the band_data_to_dict functions, which are called by this function and their outputs merged. Only the master process returns meaningful data.

Some details on the conventions for the returned data:

  • ρ: (fftsize[1], fftsize[2], fftsize[3], nspin) array of density on real-space grid.
  • energies: Dictionary / subdirectory containing the energy terms
  • converged: Has the SCF reached convergence
  • norm_Δρ: Most recent change in ρ during an SCF step
  • occupation_threshold: Threshold below which orbitals are considered unoccupied
  • n_bands_converge: Number of bands that have been fully converged numerically.
  • n_iter: Number of iterations.
source
DFTK.select_eigenpairs_all_kblocksMethod
select_eigenpairs_all_kblocks(eigres, range) -> NamedTuple

Function to select a subset of eigenpairs on each $k$-Point. Works on the Tuple returned by diagonalize_all_kblocks.

source
DFTK.self_consistent_fieldMethod
self_consistent_field(
    basis::PlaneWaveBasis{T, VT} where VT<:Real;
    ρ,
    ψ,
    tol,
    is_converged,
    maxiter,
    mixing,
    damping,
    solver,
    eigensolver,
    diagtolalg,
    nbandsalg,
    fermialg,
    callback,
    compute_consistent_energies,
    response
) -> NamedTuple{(:ham, :basis, :energies, :ρ, :eigenvalues, :occupation, :εF, :ψ, :response, :converged, :occupation_threshold, :α, :n_iter, :n_bands_converge, :diagonalization, :stage, :algorithm, :runtime_ns), <:Tuple{Hamiltonian, PlaneWaveBasis{T, VT} where {T<:ForwardDiff.Dual, VT<:Real}, Energies, Vararg{Any, 15}}}
self_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 orbitals
  • tol: Tolerance for the density change ($\|ρ_\text{out} - ρ_\text{in}\|$) to flag convergence. Default is 1e-6.
  • is_converged: Convergence control callback. Typical objects passed here are ScfConvergenceDensity(tol) (the default), ScfConvergenceEnergy(tol) or ScfConvergenceForce(tol).
  • maxiter: Maximal number of SCF iterations
  • mixing: Mixing method, which determines the preconditioner $P^{-1}$ in the above equation. Typical mixings are LdosMixing, KerkerMixing, SimpleMixing or DielectricMixing. Default is LdosMixing()
  • damping: Damping parameter $α$ in the above equation. Default is 0.8.
  • nbandsalg: By default DFTK uses nbandsalg=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 a FixedBands or AdaptiveBands. Beware that with non-zero temperature, the convergence of the SCF algorithm may be limited by the default_occupation_threshold() parameter. For highly accurate calculations we thus recommend increasing the occupation_threshold of the AdaptiveBands.
  • callback: Function called at each SCF iteration. Usually takes care of printing the intermediate state.
source
DFTK.simpsonFunction
simpson(x, y)

Integrate y(x) over x using Simpson's method quadrature.

source
DFTK.solve_ΩplusKMethod
solve_ΩplusK(
    basis::PlaneWaveBasis{T, VT} where VT<:Real,
    ψ,
    rhs,
    occupation;
    callback,
    tol
) -> NamedTuple{(:δψ, :converged, :tol, :residual_norm, :n_iter), <:Tuple{Any, Bool, Float64, Any, Int64}}
solve_ΩplusK(basis::PlaneWaveBasis{T}, ψ, res, occupation;
             tol=1e-10, verbose=false) where {T}

Return δψ where (Ω+K) δψ = rhs

source
DFTK.solve_ΩplusK_splitMethod
solve_ΩplusK_split(
    ham::Hamiltonian,
    ρ::AbstractArray{T},
    ψ,
    occupation,
    εF,
    eigenvalues,
    rhs;
    tol,
    tol_sternheimer,
    verbose,
    occupation_threshold,
    q,
    kwargs...
)

Solve 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ψ).

source
DFTK.spglib_standardize_cellMethod
spglib_standardize_cell(
    lattice::AbstractArray{T},
    atom_groups,
    positions;
    ...
) -> NamedTuple{(:lattice, :atom_groups, :positions, :magnetic_moments), <:Tuple{Matrix, Any, Any, Vector{StaticArraysCore.SVector{3, Float64}}}}
spglib_standardize_cell(
    lattice::AbstractArray{T},
    atom_groups,
    positions,
    magnetic_moments;
    correct_symmetry,
    primitive,
    tol_symmetry
) -> NamedTuple{(:lattice, :atom_groups, :positions, :magnetic_moments), <:Tuple{Matrix, Any, Any, Vector{StaticArraysCore.SVector{3, Float64}}}}

Returns 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.

source
DFTK.sphericalbesselj_fastMethod
sphericalbesselj_fast(l::Integer, x) -> Any
sphericalbesselj_fast(l::Integer, x::Number)

Returns the spherical Bessel function of the first kind $j_l(x)$. Consistent with Wikipedia and with SpecialFunctions.sphericalbesselj. Specialized for integer $0 ≤ l ≤ 5$.

source
DFTK.spin_componentsMethod
spin_components(
    spin_polarization::Symbol
) -> Union{Bool, Tuple{Symbol}, Tuple{Symbol, Symbol}}

Explicit spin components of the KS orbitals and the density

source
DFTK.split_evenlyMethod
split_evenly(itr, N) -> Any

Split an iterable evenly into N chunks, which will be returned.

source
DFTK.standardize_atomsFunction
standardize_atoms(
    lattice,
    atoms,
    positions;
    ...
) -> NamedTuple{(:lattice, :atoms, :positions, :magnetic_moments), <:Tuple{Matrix, Any, Any, Vector{StaticArraysCore.SVector{3, Float64}}}}
standardize_atoms(
    lattice,
    atoms,
    positions,
    magnetic_moments;
    kwargs...
) -> NamedTuple{(:lattice, :atoms, :positions, :magnetic_moments), <:Tuple{Matrix, Any, Any, Vector{StaticArraysCore.SVector{3, Float64}}}}

Apply 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.

source
DFTK.symmetries_preserving_kgridMethod
symmetries_preserving_kgrid(symmetries, kcoords) -> Any

Filter out the symmetry operations that don't respect the symmetries of the discrete BZ grid

source
DFTK.symmetries_preserving_rgridMethod
symmetries_preserving_rgrid(symmetries, fft_size) -> Any

Filter out the symmetry operations that don't respect the symmetries of the discrete real-space grid

source
DFTK.symmetrize_forcesMethod
symmetrize_forces(model::Model, forces; symmetries)

Symmetrize the forces in reduced coordinates, forces given as an array forces[iel][α,i].

source
DFTK.symmetrize_stressesMethod
symmetrize_stresses(model::Model, stresses; symmetries)

Symmetrize the stress tensor, given as a Matrix in Cartesian coordinates

source
DFTK.symmetrize_ρMethod
symmetrize_ρ(
    basis,
    ρ::AbstractArray{T};
    symmetries,
    do_lowpass
) -> Any

Symmetrize a density by applying all the basis (by default) symmetries and forming the average.

source
DFTK.symmetry_operationsFunction
symmetry_operations(
    lattice,
    atoms,
    positions;
    ...
) -> Union{Vector{SymOp{Bool}}, Vector{SymOp{Float64}}}
symmetry_operations(
    lattice,
    atoms,
    positions,
    magnetic_moments;
    tol_symmetry,
    check_symmetry
) -> Union{Vector{SymOp{Bool}}, Vector{SymOp{Float64}}}

Return the symmetries given an atomic structure with optionally designated magnetic moments on each of the atoms. The symmetries are determined using spglib.

source
DFTK.symmetry_operationsMethod
symmetry_operations(
    hall_number::Integer
) -> Vector{SymOp{Float64}}

Return 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.

source
DFTK.synchronize_deviceMethod
synchronize_device(_::DFTK.AbstractArchitecture)

Synchronize 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).

source
DFTK.to_cpuMethod
to_cpu(x::AbstractArray) -> Array

Transfer an array from a device (typically a GPU) to the CPU.

source
DFTK.to_deviceMethod
to_device(_::DFTK.CPU, x) -> Any

Transfer an array to a particular device (typically a GPU)

source
DFTK.todictMethod
todict(energies::Energies) -> Dict

Convert an Energies struct to a dictionary representation

source
DFTK.todictMethod
todict(model::Model) -> Dict{String, Any}

Convert a Model struct to a dictionary representation. Intended to give a condensed set of useful metadata to post-processing scripts or for storing computational results (e.g. bands, bloch waves etc.).

Some details on the conventions for the returned data:

  • lattice, recip_lattice: Always a zero-padded 3x3 matrix, independent on the actual dimension
  • atomic_positions, atomic_positions_cart: Atom positions in fractional or Cartesian coordinates, respectively.
  • atomic_symbols: Atomic symbols if known.
  • terms: Some rough information on the terms used for the computation.
  • n_electrons: Number of electrons, may be missing if εF is fixed instead
  • εF: Fixed Fermi level to use, may be missing if n_electrons is is specified instead.
source
DFTK.todictMethod
todict(basis::PlaneWaveBasis) -> Dict{String, Any}

Convert a PlaneWaveBasis struct to a dictionary representation. Intended to give a condensed set of useful metadata to post-processing scripts or for storing computational results (e.g. bands, bloch waves etc.). As such the function is lossy and might not keep all data consistently. Returns the same result on all MPI processors. See also the todict function for the Model, which is called from this one to merge the data of both outputs.

Some details on the conventions for the returned data:

  • dvol: Volume element for real-space integration
  • variational: Is the k-point specific basis (for ψ) variationally consistent with the basis for ρ.
  • kweights: Weights for the k-points, summing to 1.0
source
DFTK.total_local_potentialMethod
total_local_potential(ham::Hamiltonian) -> Any

Get the total local potential of the given Hamiltonian, in real space in the spin components.

source
DFTK.transfer_blochwaveMethod
transfer_blochwave(
    ψ_in,
    basis_in::PlaneWaveBasis,
    basis_out::PlaneWaveBasis
) -> Vector

Transfer Bloch wave between two basis sets. Limited feature set.

source
DFTK.transfer_blochwave_kptMethod
transfer_blochwave_kpt(
    ψk_in,
    basis::PlaneWaveBasis,
    kpt_in,
    kpt_out,
    ΔG
) -> Any

Transfer 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.

source
DFTK.transfer_blochwave_kptMethod
transfer_blochwave_kpt(
    ψk_in,
    basis_in::PlaneWaveBasis,
    kpt_in::Kpoint,
    basis_out::PlaneWaveBasis,
    kpt_out::Kpoint
) -> Any

Transfer an array ψk defined on basisin $k$-point kptin to basisout $k$-point kptout.

source
DFTK.transfer_densityMethod
transfer_density(
    ρ_in,
    basis_in::PlaneWaveBasis{T, VT} where VT<:Real,
    basis_out::PlaneWaveBasis{T, VT} where VT<:Real
) -> Any

Transfer 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.

source
DFTK.transfer_mappingMethod
transfer_mapping(
    basis_in::PlaneWaveBasis,
    kpt_in::Kpoint,
    basis_out::PlaneWaveBasis,
    kpt_out::Kpoint
) -> Tuple{Any, Any}

Compute 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).

source
DFTK.transfer_mappingMethod
transfer_mapping(
    basis_in::PlaneWaveBasis,
    basis_out::PlaneWaveBasis
) -> Base.Iterators.Zip{Tuple{Array{CartesianIndices{3, Tuple{UnitRange{Int64}, UnitRange{Int64}, UnitRange{Int64}}}, 3}, Array{CartesianIndices{3, Tuple{UnitRange{Int64}, UnitRange{Int64}, UnitRange{Int64}}}, 3}}}

Compute 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).

source
DFTK.trapezoidalFunction
trapezoidal(x, y)

Integrate y(x) over x using trapezoidal method quadrature.

source
DFTK.unfold_bzMethod
unfold_bz(basis::PlaneWaveBasis) -> PlaneWaveBasis

" 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).

source
DFTK.versioninfoFunction
versioninfo()
versioninfo(io::IO)
DFTK.versioninfo([io::IO=stdout])

Summary of version and configuration of DFTK and its key dependencies.

source
DFTK.weighted_ksumMethod
weighted_ksum(basis::PlaneWaveBasis, array) -> Any

Sum an array over kpoints, taking weights into account

source
DFTK.write_w90_eigMethod
write_w90_eig(fileprefix::String, eigenvalues; n_bands)

Write the eigenvalues in a format readable by Wannier90.

source
DFTK.write_w90_winMethod
write_w90_win(
    fileprefix::String,
    basis::PlaneWaveBasis;
    bands_plot,
    wannier_plot,
    kwargs...
)

Write a win file at the indicated prefix. Parameters to Wannier90 can be added as kwargs: e.g. num_iter=500.

source
DFTK.write_wannier90_filesMethod
write_wannier90_files(
    preprocess_call,
    scfres;
    n_bands,
    n_wannier,
    projections,
    fileprefix,
    wannier_plot,
    kwargs...
)

Shared file writing code for Wannier.jl and Wannier90.

source
DFTK.ylm_realMethod
ylm_real(
    l::Integer,
    m::Integer,
    rvec::AbstractArray{T, 1}
) -> Any

Returns the $(l,m)$ real spherical harmonic $Y_l^m(r)$. Consistent with Wikipedia.

source
DFTK.zeros_likeFunction
zeros_like(X::AbstractArray) -> Any
zeros_like(
    X::AbstractArray,
    T::Type,
    dims::Integer...
) -> Any

Create 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.

source
DFTK.@timingMacro

Shortened version of the @timeit macro from TimerOutputs, which writes to the DFTK timer.

source
DFTK.Smearing.occupationFunction
occupation(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.

source