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

DFTK –- The density-functional toolkit. Provides functionality for experimenting with plane-wave density-functional theory algorithms.

source
DFTK.AtomicNonlocalType

Nonlocal term coming from norm-conserving pseudopotentials in Kleinmann-Bylander form. $\text{Energy} = \sum_a \sum_{ij} \sum_{n} f_n <ψ_n|p_{ai}> D_{ij} <p_{aj}|ψ_n>.$

source
DFTK.DielectricMixingType

We use a simplification of the Resta model DOI 10.1103/physrevb.16.2717 and set $χ_0(q) = \frac{C_0 G^2}{4π (1 - C_0 G^2 / k_{TF}^2)}$ where $C_0 = 1 - ε_r$ with $ε_r$ being the macroscopic relative permittivity. We neglect $K_\text{xc}$, such that $J^{-1} ≈ \frac{k_{TF}^2 - C_0 G^2}{ε_r k_{TF}^2 - C_0 G^2}$

By default it assumes a relative permittivity of 10 (similar to Silicon). εr == 1 is equal to SimpleMixing and εr == Inf to KerkerMixing. The mixing is applied to $ρ$ and $ρ_\text{spin}$ in the same way.

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

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

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

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

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

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 = sum{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.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.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 arXiv:2009.01665

source
DFTK.MagneticType

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

source
DFTK.ModelMethod
Model(lattice; n_electrons, atoms, magnetic_moments, terms, temperature,
               smearing, spin_polarization, symmetry)

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

G_to_r and r_to_G convert between these representations.

source
DFTK.PlaneWaveBasisMethod

Creates a PlaneWaveBasis using the kinetic energy cutoff Ecut and a Monkhorst-Pack $k$-point grid. The MP grid can either be specified directly with kgrid providing the number of points in each dimension and kshift the shift (0 or 1/2 in each direction). If not specified a grid is generated using kgrid_from_minimal_spacing with a minimal spacing of 2π * 0.022 per Bohr.

If use_symmetry is true (default) the symmetries of the crystal are used to reduce the number of $k$-points which are treated explicitly. In this case all guess densities and potential functions must agree with the crystal symmetries or the result is undefined.

source
DFTK.PreconditionerTPAType

(simplified version of) Tetter-Payne-Allan preconditioning ↑ M.P. Teter, M.C. Payne and D.C. Allan, Phys. Rev. B 40, 12255 (1989).

source
DFTK.PspHghMethod
PspHgh(Zion::Number, rloc::Number, cloc::Vector, rp::Vector, h::Vector;
       identifier="", description="")

Construct a Hartwigsen, Goedecker, Teter, Hutter separable dual-space Gaussian pseudopotential (1998). The required parameters are the ionic charge Zion (total charge - valence electrons), the range for the local Gaussian charge distribution rloc, the coefficients for the local part cloc, the projector radius rp (one per AM channel) and the non-local coupling coefficients between the projectors h (one matrix per AM channel).

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.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
DFTK.CROPFunction

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_to_rMethod
G_to_r(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
DFTK.G_vectorsMethod

Return the list of wave vectors (integer coordinates) for the cubic basis set.

source
DFTK.G_vectorsMethod

The list of G vectors of a given basis or kpoint, in reduced coordinates.

source
DFTK.HybridMixingMethod

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 arXiv:2009.01665

Important kwargs passed on to χ0Mixing

  • RPA: Is the random-phase approximation used for the kernel (i.e. only Hartree kernel is used and not XC kernel)
  • verbose: Run the GMRES in verbose mode.
source
DFTK.IncreaseMixingTemperatureMethod

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

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

Determine the tolerance used for the next diagonalization. This function takes $|ρnext - ρin|$ and multiplies it with ratio_ρdiff to get the next diagtol, ensuring additionally that the returned value is between diagtol_min and diagtol_max and never increases.

source
DFTK.ScfPlotTraceFunction

Plot the trace of an SCF, i.e. the absolute error of the total energy at each iteration versus the converged energy in a semilog plot. By default a new plot canvas is generated, but an existing one can be passed and reused along with kwargs for the call to plot!. Requires Plots to be loaded.

source
DFTK.apply_KMethod
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, δρ; kwargs...)

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

source
DFTK.apply_ksymopMethod

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_ΩMethod
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

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

Note: This function assumes that all bands contained in ψ and eigenvalues are sufficiently converged. By default the self_consistent_field routine of DFTK returns 3 extra bands, which are not converged by the eigensolver (see n_ep_extra parameter). These should be discarded before using this function.

source
DFTK.ase_atoms_translation_mapMethod

Map translating between the atoms datastructure used in DFTK and the indices of the respective atoms in the ASE.Atoms object in pyobj.

source
DFTK.atom_decay_lengthMethod

Get the lengthscale of the valence density for an atom with n_elec_core core and n_elec_valence valence electrons.

source
DFTK.build_fft_plansMethod

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_projection_vectors_Method

Build projection vectors for a atoms array generated by term_nonlocal

Hat = sumij Cij |pi> <pj| Hper = sumR sumij Cij |pi(x-R)> <pj(x-R)| = sumR sum_ij Cij |pi(x-R)> <pj(x-R)|

<ekG'|Hper|ekG> = ... = 1/Ω sumij Cij pihat(k+G') pjhat(k+G)^*

where pihat(q) = ∫_R^3 pi(r) e^{-iqr} dr

We store 1/√Ω pihat(k+G) in proj_vectors.

source
DFTK.bzmesh_ir_wedgeMethod
 bzmesh_ir_wedge(kgrid_size, symmetries; kshift=[0, 0, 0])

Construct the irreducible wedge of a uniform Brillouin zone mesh for sampling $k$-points. The function returns a tuple (kcoords, ksymops), where kcoords are the list of irreducible $k$-points and ksymops are a list of symmetry operations for regenerating the full mesh. symmetries is the tuple returned from symmetry_operations(lattice, atoms, magnetic_moments). tol_symmetry is the tolerance used for searching for symmetry operations.

source
DFTK.bzmesh_uniformMethod
bzmesh_uniform(kgrid_size; kshift=[0, 0, 0])

Construct a (shifted) uniform Brillouin zone mesh for sampling the $k$-points. The function returns a tuple (kcoords, ksymops), where kcoords are the list of $k$-points and ksymops are a list of symmetry operations (for interface compatibility with PlaneWaveBasis and bzmesh_irreducible. No symmetry reduction is attempted, such that there will be prod(kgrid_size) $k$-points returned and all symmetry operations are the identity.

source
DFTK.compute_Ak_gaussian_guessMethod

Compute the matrix $[A_k]_{m,n} = \langle ψ_m^k | g^{\text{per}}_n \rangle$

$g^{per}_n$ are periodized gaussians whose respective centers are given as an (num_bands,1) array [ [center 1], ... ].

Centers are to be given in lattice coordinates and G_vectors in reduced coordinates. The dot product is computed in the Fourier space.

Given an orbital $g_n$, the periodized orbital is defined by : $g^{per}_n = \sum\limits_{R \in {\rm lattice}} g_n( \cdot - R)$. The Fourier coefficient of $g^{per}_n$ at any G is given by the value of the Fourier transform of $g_n$ in G.

source
DFTK.compute_densityMethod
compute_density(basis::PlaneWaveBasis, ψ::AbstractVector, occupation::AbstractVector)

Compute the density and spin 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. If the Model underlying the basis is not collinear the spin density is nothing.

source
DFTK.compute_fft_sizeMethod

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.

source
DFTK.compute_forcesMethod

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 [[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_forces_cartMethod

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_kernelMethod
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_nosMethod
compute_nos(ε, basis, eigenvalues; smearing=basis.model.smearing,
            temperature=basis.model.temperature)

The number of Kohn-Sham states in a temperature window of width temperature around the energy ε contributing to the DOS at temperature T.

This quantity is not a physical quantity, but rather a dimensionless approximate measure for how well properties near the Fermi surface are sampled with the passed smearing and temperature T. It increases with both T and better sampling of the BZ with $k$-points. A value $\gg 1$ indicates a good sampling of properties near the Fermi surface.

source
DFTK.compute_occupation_bandgapMethod

Find Fermi level and occupation for the given parameters, assuming a band gap and zero temperature. This function is for DEBUG purposes only, and the finite-temperature version with 0 temperature should be preferred.

source
DFTK.compute_recip_latticeMethod

Compute the reciprocal lattice. Takes special care of 1D or 2D cases. 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_transfer_matrixMethod

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_χ0Method

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

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

Computes the ground state by direct minimization. kwargs... are passed to Optim.Options(). Note that the resulting ψ are not necessarily eigenvectors of the Hamiltonian.

source
DFTK.disable_threadingMethod

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

source
DFTK.divergence_realMethod

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_ewaldMethod

Compute the electrostatic interaction energy per unit cell between point charges in a uniform background of compensating charge to yield net neutrality. the lattice and recip_lattice should contain the lattice and reciprocal lattice vectors as columns. charges and positions are the point charges and their positions (as an array of arrays) in fractional coordinates. If forces is not nothing, minus the derivatives of the energy with respect to positions is computed.

source
DFTK.energy_psp_correctionMethod
energy_psp_correction(model)

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

source
DFTK.eval_psp_local_fourierMethod
eval_psp_local_fourier(psp, q)

Evaluate the local part of the pseudopotential in reciprocal space: V(q) = ∫R^3 Vloc(r) e^{-iqr} dr = 4π ∫{R+} sin(qr)/q r e^{-iqr} dr

source
DFTK.eval_psp_projector_fourierMethod
eval_psp_projector_fourier(psp, i, l, q)

Evaluate the radial part of the i-th projector for angular momentum l at the reciprocal vector with modulus q: p(q) = ∫R^3 p{il}(r) e^{-iqr} dr = 4π ∫{R+} r^2 p{il}(r) j_l(q r) dr

source
DFTK.eval_psp_projector_realMethod
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.gather_kptsMethod

Gather the distributed data of a quantity depending on k-Points on the master process and return it. On the other (non-master) processes nothing is returned.

source
DFTK.gather_kptsMethod

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 to extract data for post-processing and serialisation to disk.

source
DFTK.gaussian_superpositionMethod

Build a superposition of Gaussians as a guess for the density and magnetisation. Expects a list of tuples (coefficient, length, position) for each of the Gaussian, which follow the functional form

\[\hat{ρ}(G) = \text{coefficient} \exp\left(-(2π \text{length} |G|)^2\right)\]

and are placed at position (in fractional coordinates).

source
DFTK.get_spglib_latticeMethod

Returns crystallographic conventional cell according to the International Table of Crystallography Vol A (ITA) in case to_primitive=false. If to_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.guess_densityFunction
guess_density(basis, magnetic_moments)

Build a superposition of atomic densities (SAD) guess density.

We take for the guess density a gaussian centered around the atom, of length specified by atom_decay_length, normalized to get the right number of electrons

\[\hat{ρ}(G) = Z \exp\left(-(2π \text{length} |G|)^2\right)\]

When magnetic moments are provided, construct a symmetry-broken density guess. The magnetic moments should be specified in units of $μ_B$.

source
DFTK.high_symmetry_kpathMethod

Extract the high-symmetry $k$-point path corresponding to the passed model using Brillouin.jl. 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. The kline_density is given in number of $k$-points per inverse bohrs (i.e. overall in units of length).

Issues a warning in case the passed lattice does not match the expected primitive.

source
DFTK.index_G_vectorsMethod

Return the index tuple I such that G_vectors(basis)[I] == G or the index i such that G_vectors(kpoint)[i] == G. Returns nothing if outside the range of valid wave vectors.

source
DFTK.interpolate_densityMethod

Interpolate a function expressed in a basis basis_in to a basis basis_out This interpolation uses a very basic real-space algorithm, and makes a DWIM-y attempt to take into account the fact that basisout can be a supercell of basisin

source
DFTK.interpolate_kpointMethod

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

source
DFTK.is_metalFunction
is_metal(band_data, ε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.kgrid_from_minimal_n_kpointsMethod

Selects a kgrid size which ensures that at least a n_kpoints total number of $k$-points are used. The distribution of $k$-points amongst coordinate directions is as uniformly as possible, trying to achieve an identical minimal spacing in all directions.

source
DFTK.kgrid_from_minimal_spacingMethod

Selects a kgrid size to ensure a minimal spacing (in inverse Bohrs) between kpoints. A reasonable spacing is 0.13 inverse Bohrs (around $2π * 0.04 \AA^{-1}$).

source
DFTK.lda_c_vwn!Method

LDA correlation according to Vosko Wilk,and Nusair, (DOI 10.1139/p80-159)

source
DFTK.lda_x!Method

LDA Slater exchange (DOI: 10.1017/S0305004100016108 and 10.1007/BF01340281)

source
DFTK.list_pspFunction
list_psp(element; functional, family, core, datadir_psp)

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_atomsMethod

Return a DFTK-compatible atoms object loaded from an ASE Atoms, a pymatgen Structure or a compatible file (e.g. xyz file, cif file etc.)

source
DFTK.load_atoms_pymatgenMethod

Load a DFTK-compatible atoms representation from a supported pymatgen object. All atoms are using a Coulomb model.

source
DFTK.load_latticeMethod

Return a DFTK-compatible lattice object loaded from an ASE Atoms, a pymatgen Structure or a compatible file (e.g. xyz file, cif file etc.)

source
DFTK.load_latticeMethod

Load a DFTK-compatible lattice object from a supported python object (e.g. pymatgen or ASE)

source
DFTK.load_pspMethod

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

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

source
DFTK.model_LDAMethod

Build an LDA model (Teter93 parametrization) from the specified atoms.

source
DFTK.model_atomicMethod

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

source
DFTK.mpi_nprocsFunction

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

source
DFTK.newtonMethod
newton(basis::PlaneWaveBasis{T}; ψ0=nothing,
       tol=1e-6, tol_cg=1e-10, maxiter=20, verbose=false,
       callback=NewtonDefaultCallback(),
       is_converged=NewtonConvergenceDensity(tol))

Newton algorithm. Be careful that the starting point needs to be not too far from the solution.

source
DFTK.overlap_Mmn_k_kpbMethod

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

G_shift is the "shifting" vector, correction due to the periodicity conditions imposed on k -> ψk. It is non zero if kpb is taken in another unit cell of the reciprocal lattice. We use here that : ``u{n(k + Gshift)}(r) = e^{-i*\langle Gshift,r \rangle} u_{nk}``

source
DFTK.plot_bandstructureMethod

Compute and plot the band structure. n_bands selects the number of bands to compute. If this value is absent and an scfres is used to start the calculation a default of n_bands_scf + 5sqrt(n_bands_scf) is used. The unit used to plot the bands can be selected using the unit parameter. Like in the rest of DFTK Hartree is used by default. Another standard choices is unit=u"eV" (electron volts). The kline_density is given in number of $k$-points per inverse bohrs (i.e. overall in units of length).

source
DFTK.psp_local_polynomialFunction

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

source
DFTK.psp_projector_polynomialFunction

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

source
DFTK.qcut_psp_localMethod

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

source
DFTK.qcut_psp_projectorMethod

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

source
DFTK.r_to_G!Method

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

source
DFTK.r_to_GMethod
r_to_G(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
DFTK.r_vectorsMethod

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

source
DFTK.read_w90_nnkpMethod

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

Run an SCF in ABINIT starting from a DFTK Model and some extra parameters. Write the result to the output directory in ETSF Nanoquanta format and return the EtsfFolder object.

source
DFTK.run_abinit_scfMethod

Run an SCF in ABINIT starting from the input file infile represented as a abipy.abilab.AbinitInput python object. Write the result to the output directory in ETSF Nanoquanta format and return the result as an EtsfFolder object.

source
DFTK.run_wannier90Method

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 using random Gaussians as guesses. All keyword arguments supported by Wannier90 for the disentanglement may be added as keyword arguments. The function returns the fileprefix.

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_scfresMethod
save_scfres(filename, scfres)

Save an scfres obtained from self_consistent_field to a file. The format is determined from the file extension. Currently the following file extensions are recognized and supported:

  • jld2: A JLD2 file. Stores the complete state and can be used (with load_scfres) to restart an SCF from a checkpoint or post-process an SCF solution. See Saving SCF results on disk and SCF checkpoints for details.
  • vts: A VTK file for visualisation e.g. in paraview. Stores the density, spin density and some metadata (energy, Fermi level, occupation etc.). Supports these keyword arguments:
    • save_ψ: Save the real-space representation of the orbitals as well (may lead to larger files).
    • extra_data: Dict{String,Array} with additional data on the 3D real-space grid to store into the VTK file.
No compatibility guarantees

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

source
DFTK.scf_damping_quadratic_modelMethod

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_nlsolve_solverFunction

Create a NLSolve-based SCF solver, by default using an Anderson-accelerated fixed-point scheme, keeping m steps for Anderson acceleration. See the NLSolve documentation for details about the other parameters and methods.

source
DFTK.solve_ΩplusKMethod
solve_ΩplusK(basis::PlaneWaveBasis{T}, ψ, res, occupation;
             tol_cg=1e-10, verbose=false) where T

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

source
DFTK.spglib_atomsFunction

Convert the DFTK atoms datastructure into a tuple of datastructures for use with spglib. positions contains positions per atom, numbers contains the mapping atom to a unique number for each indistinguishable element, spins contains the $z$-component of the initial magnetic moment on each atom, mapping contains the mapping of the numbers to the element objects in DFTK and collinear whether the atoms mark a case of collinear spin or not. Notice that if collinear is false then spins is garbage.

source
DFTK.standardize_atomsFunction

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.symmetrize_ρMethod

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

source
DFTK.transfer_mappingMethod

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

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

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

source
DFTK.ylm_realMethod

Returns the (l,m) real spherical harmonic Ylm(r). Consistent with https://en.wikipedia.org/wiki/Tableofsphericalharmonics#Realsphericalharmonics

source
DFTK.@timingMacro

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

source
DFTK.@timing_seqMacro

Similar to @timing, but disabled in parallel runs. Should be used to time threaded regions, since TimerOutputs is not thread-safe and breaks otherwise.

source
DFTK.Smearing.entropyMethod

Entropy. Note that this is a function of the energy x, not of occupation(x). This function satisfies s' = x f' (see https://www.vasp.at/vasp-workshop/k-points.pdf p. 12 and https://arxiv.org/pdf/1805.07144.pdf p. 18.

source