API reference
This page provides a plain list of all documented functions, structs, modules and macros in DFTK. Note that this list is neither structured, complete nor particularly clean, so it only provides rough orientation at the moment. The best reference is the code itself.
DFTK.DFTK
— ModuleDFTK –- The density-functional toolkit. Provides functionality for experimenting with plane-wave density-functional theory algorithms.
DFTK.DFTK_DATADIR
— ConstantThe default search location for Pseudopotential data files
DFTK.timer
— ConstantTimerOutput object used to store DFTK timings.
DFTK.Applyχ0Model
— TypeFull χ0 application, optionally dropping terms or disabling Sternheimer. All keyword arguments passed to apply_χ0
.
DFTK.AtomicLocal
— TypeAtomic local potential defined by model.atoms
.
DFTK.AtomicNonlocal
— TypeNonlocal term coming from norm-conserving pseudopotentials in Kleinmann-Bylander form. $\text{Energy} = \sum_a \sum_{ij} \sum_{n} f_n <ψ_n|p_{ai}> D_{ij} <p_{aj}|ψ_n>.$
DFTK.DielectricMixing
— TypeWe use a simplification of the Resta model DOI 10.1103/physrevb.16.2717 and set $χ_0(q) = \frac{C_0 G^2}{4π (1 - C_0 G^2 / k_{TF}^2)}$ where $C_0 = 1 - ε_r$ with $ε_r$ being the macroscopic relative permittivity. We neglect $K_\text{xc}$, such that $J^{-1} ≈ α \frac{k_{TF}^2 - C_0 G^2}{ε_r k_{TF}^2 - C_0 G^2}$
By default it assumes a relative permittivity of 10 (similar to Silicon). εr == 1
is equal to SimpleMixing
and εr == Inf
to KerkerMixing
. The mixing is applied to $ρ$ and $ρ_\text{spin}$ in the same way.
DFTK.DielectricModel
— TypeA localised dielectric model for $χ_0$:
where $C_0 = 1 - ε_r$, L(r)
is a real-space localisation function and otherwise the same conventions are used as in DielectricMixing
.
DFTK.ElementCohenBergstresser
— MethodElement where the interaction with electrons is modelled as in CohenBergstresser1966. Only the homonuclear lattices of the diamond structure are implemented (i.e. Si, Ge, Sn).
key
may be an element symbol (like :Si
), an atomic number (e.g. 14
) or an element name (e.g. "silicon"
)
DFTK.ElementCoulomb
— MethodElement interacting with electrons via a bare Coulomb potential (for all-electron calculations) key
may be an element symbol (like :Si
), an atomic number (e.g. 14
) or an element name (e.g. "silicon"
)
DFTK.ElementPsp
— MethodElement interacting with electrons via a pseudopotential model. key
may be an element symbol (like :Si
), an atomic number (e.g. 14
) or an element name (e.g. "silicon"
)
DFTK.Energies
— TypeA simple struct to contain a vector of energies, and utilities to print them in a nice format.
DFTK.Entropy
— TypeEntropy term -TS, where S is the electronic entropy. Turns the energy E into the free energy F=E-TS. This is in particular useful because the free energy, not the energy, is minimized at self-consistency.
DFTK.EtsfFolder
— MethodInitialize a EtsfFolder from the path to the folder which contains the data in the ETSF Nanoquanta format.
DFTK.Ewald
— TypeEwald term: electrostatic energy per unit cell of the array of point charges defined by model.atoms
in a uniform background of compensating charge yielding net neutrality.
DFTK.ExternalFromFourier
— TypeExternal potential from the (unnormalized) Fourier coefficients V(G)
G is passed in cartesian coordinates
DFTK.ExternalFromReal
— TypeExternal potential from an analytic function V
(in cartesian coordinates). No low-pass filtering is performed.
DFTK.FourierMultiplication
— TypeFourier space multiplication, like a kinetic energy term: (Hψ)(G) = multiplier(G) ψ(G)
DFTK.Hartree
— TypeHartree term: for a decaying potential V the energy would be
1/2 ∫ρ(x)ρ(y)V(x-y) dxdy
with the integral on x in the unit cell and of y in the whole space. For the Coulomb potential with periodic boundary conditions, this is rather
1/2 ∫ρ(x)ρ(y) G(x-y) dx dy
where G is the Green's function of the periodic Laplacian with zero mean (-Δ G = sum{R} 4π δR, integral of G zero on a unit cell).
DFTK.KerkerDosMixing
— TypeThe same as KerkerMixing
, but the Thomas-Fermi wavevector is computed from the current density of states at the Fermi level.
DFTK.KerkerMixing
— TypeKerker mixing: $J^{-1} ≈ \frac{α |G|^2}{k_{TF}^2 + |G|^2}$ where $k_{TF}$ is the Thomas-Fermi wave vector. For spin-polarized calculations by default the spin density is not preconditioned. Unless a non-default value for $ΔDOS$ is specified. This value should roughly be the expected difference in density of states (per unit volume) between spin-up and spin-down.
Notes:
- Abinit calls $1/k_{TF}$ the dielectric screening length (parameter dielng)
DFTK.Kinetic
— TypeKinetic energy: 1/2 sumn fn ∫ |∇ψn|^2.
DFTK.Kpoint
— TypeDiscretization information for kpoint-dependent quantities such as orbitals. More generally, a kpoint is a block of the Hamiltonian; eg collinear spin is treated by doubling the number of kpoints.
DFTK.LdosModel
— TypeRepresents the LDOS-based $χ_0$ model
where $D_\text{loc}$ is the local density of states and $D$ the density of states. For details see Herbst, Levitt 2020 arXiv:2009.01665
DFTK.LibxcDensity
— TypeCompute density in real space and its derivatives starting from ρ
DFTK.Magnetic
— TypeMagnetic term $A⋅(-i∇)$. It is assumed (but not checked) that $∇⋅A = 0$.
DFTK.MagneticFieldOperator
— TypeMagnetic field operator A⋅(-i∇).
DFTK.Model
— MethodModel(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.
DFTK.NonlocalOperator
— TypeNonlocal operator in Fourier space in Kleinman-Bylander format, defined by its projectors P matrix and coupling terms D: Hψ = PDP' ψ
DFTK.NoopOperator
— TypeNoop operation: don't do anything. Useful for energy terms that don't depend on the orbitals at all (eg nuclei-nuclei interaction).
DFTK.NoopTerm
— TypeA term with a constant zero energy.
DFTK.PlaneWaveBasis
— TypeA plane-wave discretized Model
. Normalization conventions:
- Things that are expressed in the G basis are normalized so that if $x$ is the vector, then the actual function is $sum_G x_G e_G$ with $e_G(x) = e^{iG x}/sqrt(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.
DFTK.PlaneWaveBasis
— TypeCreates a new basis identical to basis
, but with a different set of kpoints
DFTK.PlaneWaveBasis
— MethodCreates a PlaneWaveBasis
using the kinetic energy cutoff Ecut
and a Monkhorst-Pack kpoint 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_size_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.
DFTK.PlaneWaveBasis
— Method" Convert a basis
into one that uses or doesn't use BZ symmetrization Mainly useful for debug purposes (e.g. in cases we don't want to bother with symmetry)
DFTK.PowerNonlinearity
— TypePower nonlinearity, with energy C ∫ρ^α where ρ is the total density
DFTK.PreconditionerNone
— TypeNo preconditioning
DFTK.PreconditionerTPA
— Type(simplified version of) Tetter-Payne-Allan preconditioning ↑ M.P. Teter, M.C. Payne and D.C. Allan, Phys. Rev. B 40, 12255 (1989).
DFTK.PspCorrection
— TypePseudopotential correction energy. TODO discuss the need for this.
DFTK.PspHgh
— MethodPspHgh(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).
DFTK.RealFourierArray
— TypeA structure to facilitate manipulations of an array of real-space type T, in both real and fourier space. Create with from_real
or from_fourier
, and access with A.real
and A.fourier
.
DFTK.RealFourierOperator
— TypeLinear operators that act on tuples (real, fourier) The main entry point is apply!(out, op, in)
which performs the operation out += op*in where out and in are named tuples (real, fourier) They also implement mul! and Matrix(op) for exploratory use.
DFTK.RealSpaceMultiplication
— TypeReal space multiplication by a potential: (Hψ)(r) V(r) ψ(r)
DFTK.SimpleMixing
— TypeSimple mixing: $J^{-1} ≈ α$
DFTK.Xc
— TypeExchange-correlation term, defined by a list of functionals and usually evaluated through libxc.
DFTK.χ0Mixing
— TypeGeneric mixing function using a model for the susceptibility composed of the sum of the χ0terms
. For valid χ0terms
See the subtypes of χ0Model
. The dielectric model is solved in real space using a GMRES. Either the full kernel (RPA=false
) or only the Hartree kernel (RPA=true
) are employed. verbose=true
lets the GMRES run in verbose mode (useful for debugging).
DFTK.CROP
— FunctionCROP-accelerated root-finding iteration for f
, starting from x0
and keeping a history of m
steps. Optionally warming
specifies the number of non-accelerated steps to perform for warming up the history.
DFTK.DOS
— MethodTotal density of states at energy ε
DFTK.G_to_r!
— MethodIn-place version of G_to_r
.
DFTK.G_to_r
— MethodG_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.
DFTK.G_vectors
— MethodReturn the list of wave vectors (integer coordinates) for the cubic basis set.
DFTK.G_vectors
— MethodThe list of G vectors of a given basis
or kpoint
, in reduced coordinates.
DFTK.G_vectors_cart
— MethodThe list of G vectors of a given basis
or kpoint
, in cartesian coordinates.
DFTK.HybridMixing
— MethodThe model for the susceptibility is
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 localisation function L(r)
. For details see Herbst, Levitt 2020 arXiv:2009.01665
Important kwargs
passed on to χ0Mixing
α
: Damping parameterRPA
: 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.
DFTK.LDOS
— MethodLocal density of states, in real space
DFTK.NOS
— MethodNOS(ε, 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.
DFTK.ScfConvergenceDensity
— MethodFlag convergence by using the L2Norm of the change between input density and unpreconditioned output density (ρout)
DFTK.ScfConvergenceEnergy
— MethodFlag convergence as soon as total energy change drops below tolerance
DFTK.ScfDefaultCallback
— MethodDefault callback function for self_consistent_field
, which prints a convergence table
DFTK.ScfDiagtol
— MethodDetermine the tolerance used for the next diagonalization. This function takes $|ρnext - ρin|$ and multiplies it with ratio_ρdiff
to get the next diagtol
, ensuring additionally that the returned value is between diagtol_min
and diagtol_max
and never increases.
DFTK.ScfPlotTrace
— FunctionPlot the trace of an SCF, i.e. the absolute error of the total energy at each iteration versus the converged energy in a semilog plot. By default a new plot canvas is generated, but an existing one can be passed and reused along with kwargs
for the call to plot!
.
DFTK.ScfSaveCheckpoints
— FunctionAdds simplistic checkpointing to a DFTK self-consistent field calculation.
DFTK.add_response_from_band!
— MethodAdds the term (f'ₙ δεₙ |ψₙ|² + 2Re fₙ ψₙ * δψₙ
to δρ_{k}
where δψₙ
is computed from δV
partly using the known, computed states and partly by solving the Sternheimer equation (if sternheimer_contribution=true
).
DFTK.anderson
— FunctionAnderson-accelerated root-finding iteration for finding a root of f
, starting from x0
and keeping a history of m
steps. Optionally warming
specifies the number of non-accelerated steps to perform for warming up the history.
DFTK.apply_kernel
— Functionapply_kernel(basis::PlaneWaveBasis, dρ, dρspin=nothing; kwargs...)
Computes the potential response to a perturbation (dρ, dρspin)
in real space. Returns the array [dV_α, dV_β]
for collinear spin-polarized calculations, else the array [dV_{tot}].
DFTK.apply_ksymop
— MethodApply a symmetry operation to eigenvectors ψk
at a given kpoint
to obtain an equivalent point in [-0.5, 0.5)^3 and associated eigenvectors (expressed in the basis of the new kpoint).
DFTK.apply_ksymop
— MethodApply a k
-point symmetry operation (the tuple (S, τ)) to a partial density.
DFTK.apply_χ0
— FunctionReturns the change in density δρ
for a given δV
. Drop all non-diagonal terms with (f(εn)-f(εm))/(εn-εm) factor less than droptol
. If sternheimer_contribution
is false, only compute excitations inside the provided orbitals.
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.
DFTK.atom_decay_length
— MethodGet the lengthscale of the valence density for an atom with n_elec_core
core and n_elec_valence
valence electrons. ```
DFTK.build_fft_plans
— MethodPlan a FFT of type T
and size fft_size
, spending some time on finding an optimal algorithm. (Inplace, out-of-place) x (forward, backward) FFT plans are returned.
DFTK.build_form_factors
— MethodBuild form factors (Fourier transforms of projectors) for an atom centered at 0.
DFTK.build_projection_vectors_
— MethodBuild projection vectors for a atoms array generated by term_nonlocal
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.
DFTK.bzmesh_ir_wedge
— Method bzmesh_ir_wedge(kgrid_size, symmetries; kshift=[0, 0, 0])
Construct the irreducible wedge of a uniform Brillouin zone mesh for sampling $k$-Points. 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.
DFTK.bzmesh_uniform
— Methodbzmesh_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.
DFTK.charge_ionic
— MethodReturn the total ionic charge of an atom type (nuclear charge - core electrons)
DFTK.charge_nuclear
— MethodReturn the total nuclear charge of an atom type
DFTK.clear_without_conjugate!
— MethodZero all elements of a Fourier space array which have no complex-conjugate partner and may thus lead to an imaginary component in real space (after an iFFT).
DFTK.compute_current
— MethodComputes the probability (not charge) current, ∑ fn Im(ψn* ∇ψn)
DFTK.compute_density
— Methodcompute_density(basis::PlaneWaveBasis, ψ::AbstractVector, occupation::AbstractVector)
Compute the density 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
.
DFTK.compute_kernel
— Methodcompute_kernel(basis::PlaneWaveBasis; kwargs...)
Computes a matrix representation of the full response kernel (derivative of potential with respect to density) in real space. For non-spin-polarized calculations the matrix dimension is prod(basis.fft_size)
× prod(basis.fft_size)
and for collinear spin-polarized cases it is 2prod(basis.fft_size)
× 2prod(basis.fft_size)
. In this case the matrix has effectively 4 blocks, which are:
i.e. corresponding to a mapping $(ρ_\text{tot}, ρ_\text{spin})^T = (ρ_α + ρ_β, ρ_α - ρ_β)^T ↦ (V_α, V_β)^T$.
DFTK.compute_partial_density!
— MethodCompute the partial density at the indicated $k$-Point and return it (in Fourier space).
DFTK.compute_χ0
— MethodCompute the independent-particle susceptibility. Will blow up for large systems. Drop all non-diagonal terms with (f(εn)-f(εm))/(εn-εm) factor less than droptol
. 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:
i.e. corresponding to a mapping ``(Vα, Vβ)^T ↦ (ρ\text{tot}, ρ\text{spin})^T = (ρα + ρβ, ρα - ρβ)^T.
DFTK.datadir_psp
— MethodReturn the data directory with pseudopotential files
DFTK.default_spin_polarization
— Method:none if no element has a magnetic moment, else :collinear or :full
DFTK.default_symmetries
— MethodDefault logic to determine the symmetry operations to be used in the model.
DFTK.determine_fft_size
— MethodDetermine the minimal grid size for the cubic basis set to be able to represent product of orbitals (with the default supersampling=2
).
Optionally optimize the grid afterwards for the FFT procedure by ensuring factorization into small primes.
The function will determine the smallest parallelepiped containing the wave vectors $|G|^2/2 \leq E_\text{cut} ⋅ \text{supersampling}^2$. For an exact representation of the density resulting from wave functions represented in the spherical basis sets, supersampling
should be at least 2
.
DFTK.diagonalize_all_kblocks
— MethodFunction for diagonalising each $k$-Point blow of ham one step at a time. Some logic for interpolating between $k$-Points is used if interpolate_kpoints
is true and if no guesses are given. eigensolver
is the iterative eigensolver that really does the work, operating on a single $k$-Block. eigensolver
should support the API eigensolver(A, X0; prec, tol, maxiter)
prec_type
should be a function that returns a preconditioner when called as prec(ham, kpt)
DFTK.direct_minimization
— MethodComputes the ground state by direct minimization. kwargs...
are passed to Optim.Options()
. Note that the resulting ψ are not necessarily eigenvectors of the Hamiltonian.
DFTK.divergence_real
— MethodCompute divergence of an operand function, which returns the cartesian x,y,z components in real space when called with the arguments 1 to 3. The divergence is also returned as a real-space array.
DFTK.energy_ewald
— MethodCompute 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.
DFTK.energy_psp_correction
— Methodenergy_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.
DFTK.eval_psp_energy_correction
— Methodeval_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.
DFTK.eval_psp_local_fourier
— Methodeval_psp_local_fourier(psp, q)
Evaluate the local part of the pseudopotential in reciprocal space.
This function computes V(q) = ∫R^3 Vloc(r) e^{-iqr} dr = 4π ∫{R+} sin(qr)/q r e^{-iqr} dr
GTH98 except they do it with plane waves normalized by 1/sqrt(Ω).
DFTK.eval_psp_local_real
— MethodDFTK.eval_psp_projection_radial
— Methodeval_psp_projection_radial(psp::PspHgh, i, l, q::Number)
Evaluate the radial part of the i
-th projector for angular momentum l
at the reciprocal vector with modulus q
.
p(q) = ∫{R+} r^2 p(r) jl(q r) dr
HGH98 except they do it with plane waves normalized by 1/sqrt(Ω).
DFTK.eval_psp_projection_radial_real
— Methodeval_psp_projection_radial_real(psp::PspHgh, i, l, q::Real)
Evaluate the radial part of the i
-th projector for angular momentum l
in real-space at the vector with modulus r
. HGH98
DFTK.fermi_level
— MethodFind the Fermi level.
DFTK.filled_occupation
— MethodMaximal occupation of a state (2 for non-spin-polarized electrons, 1 otherwise).
DFTK.find_irreducible_kpoints
— MethodImplements a primitive search to find an irreducible subset of kpoints amongst the provided kpoints.
DFTK.find_occupation
— MethodFind the occupation and Fermi level.
DFTK.find_occupation_bandgap
— MethodFind 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.
DFTK.gaussian_superposition
— MethodBuild 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
and are placed at position
(in fractional coordinates).
DFTK.guess_density
— Methodguess_density(basis)
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
DFTK.guess_spin_density
— Functionguess_spin_density(basis, magnetic_moments)
Form a spin density guess. The magnetic moments should be specified in units of $μ_B$.
DFTK.hamiltonian_with_total_potential
— MethodReturns a new Hamiltonian with local potential replaced by the given one
DFTK.index_G_vectors
— MethodReturn 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.
DFTK.interpolate_blochwave
— MethodInterpolate Bloch wave between two basis sets. Limited feature set. Currently only interpolation to a bigger grid (larger Ecut) on the same lattice supported.
DFTK.interpolate_density
— MethodInterpolate a function expressed in a basis b_in
to a basis b_out
This interpolation uses a very basic real-space algorithm, and makes a DWIM-y attempt to take into account the fact that bout can be a supercell of bin
DFTK.interpolate_kpoint
— MethodInterpolate 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
DFTK.is_metal
— Functionis_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.
DFTK.kgrid_monkhorst_pack
— MethodConstruct the coordinates of the kpoints in a (shifted) Monkorst-Pack grid
DFTK.kgrid_size_from_minimal_spacing
— FunctionSelects a kgrid_size to ensure a minimal spacing (in inverse Bohrs) between kpoints. Default is $2π * 0.04 \AA^{-1}$.
DFTK.krange_spin
— MethodReturn the index range of $k$-points that have a particular spin component.
DFTK.lda_c_vwn!
— MethodLDA correlation according to Vosko Wilk,and Nusair, (DOI 10.1139/p80-159)
DFTK.lda_x!
— MethodLDA Slater exchange (DOI: 10.1017/S0305004100016108 and 10.1007/BF01340281)
DFTK.libxc_symmetric_index
— MethodTranslation from a double index (i, j) of two symmetric axes to the linear index over the non-redundant components used in libxc
DFTK.list_psp
— Functionlist_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.
DFTK.load_atoms
— MethodLoad a DFTK-compatible atoms object from the ETSF folder. Use the scalar type T
to represent the data.
DFTK.load_atoms_pymatgen
— MethodLoad a DFTK-compatible atoms representation from a supported pymatgen object. All atoms are using a Coulomb model.
DFTK.load_basis
— MethodLoad a DFTK-compatible basis object from the ETSF folder. Use the scalar type T
to represent the data.
DFTK.load_density
— MethodLoad a DFTK-compatible density object from the ETSF folder. Use the scalar type T
to represent the data.
DFTK.load_lattice
— MethodLoad a DFTK-compatible lattice object from the ETSF folder
DFTK.load_lattice
— MethodLoad a DFTK-compatible lattice object from a supported python object (e.g. pymatgen or ASE)
DFTK.load_model
— MethodLoad a DFTK-compatible model object from the ETSF folder. Use the scalar type T
to represent the data.
DFTK.load_psp
— Methodload_psp(key; datadir_psp)
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.
DFTK.load_scfres
— Functionload_scfres(filename)
Load back an scfres
, which has previously been stored with save_scfres
. Note the warning in save_scfres
.
DFTK.local_potential_fourier
— MethodRadial local potential, in Fourier space: V(q) = int_{R^3} V(x) e^{-iqx} dx.
DFTK.local_potential_real
— MethodRadial local potential, in real space.
DFTK.model_DFT
— MethodBuild a DFT model from the specified atoms, with the specified functionals.
DFTK.model_LDA
— MethodBuild an LDA model (Teter93 parametrization) from the specified atoms.
DFTK.model_PBE
— MethodBuild an PBE-GGA model from the specified atoms.
DFTK.model_atomic
— MethodConvenience constructor, which builds a standard atomic (kinetic + atomic potential) model. Use extra_terms
to add additional terms.
DFTK.mpi_ensure_initialized
— MethodInitialize MPI. Must be called before doing any non-trivial MPI work (even in the single-process case). Unlike the MPI.Init() function, this can be called multiple times.
DFTK.mpi_nprocs
— FunctionNumber of processors used in MPI. Can be called without ensuring initialization.
DFTK.mpi_split_iterator
— MethodSplits an iterator evenly between the processes of comm
and returns the part handled by the current process.
DFTK.n_elec_core
— MethodReturn the number of core electrons
DFTK.n_elec_valence
— MethodReturn the number of valence electrons
DFTK.next_density
— MethodObtain new density ρ by diagonalizing ham
.
DFTK.normalize_kpoint_coordinate
— MethodBring kpoint coordinates into the range [-0.5, 0.5)
DFTK.parse_hgh_file
— Methodparse_hgh_file(path; identifier="")
Parse an HGH pseudopotential file and construct the PspHgh object. If identifier
is given, this identifier will be set.
DFTK.plot_bandstructure
— MethodCompute 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. Unlike the rest of DFTK bands energies are plotted in :eV
unless a different unit
is selected.
DFTK.plot_dos
— FunctionPlot the density of states over a reasonable range
DFTK.psp_local_polynomial
— FunctionThe local potential of a HGH pseudopotentials in reciprocal space can be brought to the form $Q(t) / (t^2 exp(t^2 / 2))$ where $t = r_\text{loc} q$ and Q
is a polynomial of at most degree 8. This function returns Q
.
DFTK.psp_projection_radial_polynomial
— FunctionThe nonlocal projectors of a HGH pseudopotentials in reciprocal space can be brought to the form $Q(t) exp(-t^2 / 2)$ where $t = r_l q$ and Q
is a polynomial. This function returns Q
.
DFTK.qcut_psp_local
— MethodEstimate an upper bound for the argument q
after which abs(eval_psp_local_fourier(psp, q))
is a strictly decreasing function.
DFTK.qcut_psp_projection_radial
— MethodEstimate an upper bound for the argument q
after which eval_psp_projection_radial(psp, q)
is a strictly decreasing function.
DFTK.r_to_G!
— MethodIn-place version of r_to_G!
. NOTE: If kpt
is given, not only f_fourier
but also f_real
is overwritten.
DFTK.r_to_G
— Methodr_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.
DFTK.r_vectors
— MethodReturn the list of r vectors, in reduced coordinates. By convention, this is in [0,1)^3.
DFTK.r_vectors_cart
— MethodReturn the list of r vectors, in cartesian coordinates.
DFTK.run_abinit_scf
— MethodRun 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.
DFTK.run_abinit_scf
— MethodRun 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.
DFTK.save_scfres
— Functionsave_scfres(filename, scfres)
Save an scfres
obtained from self_consistent_field
to a JLD2 file.
No guarantees are made with respect to this function at this point. It may change incompatibly between DFTK versions or stop working / be removed in the future.
DFTK.scf_damping_solver
— FunctionCreate a damped SCF solver updating the density as x = β * x_new + (1 - β) * x
DFTK.scf_nlsolve_solver
— FunctionCreate 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.
DFTK.select_eigenpairs_all_kblocks
— MethodFunction to select a subset of eigenpairs on each $k$-Point. Works on the Tuple returned by diagonalize_all_kblocks
.
DFTK.self_consistent_field
— MethodSolve the Kohn-Sham equations with a SCF algorithm, starting at ρ.
DFTK.spglib_atoms
— FunctionConvert 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.
DFTK.spin_components
— MethodExplicit spin components of the KS orbitals and the density
DFTK.split_evenly
— MethodSplit an iterable evenly into N chunks, which will be returned.
DFTK.standardize_atoms
— FunctionApply various standardisations to a lattice and a list of atoms. It uses spglib to detect symmetries (within tol_symmetry
), then cleans up the lattice according to the symmetries (unless correct_symmetry
is false
) and returns the resulting standard lattice and atoms. If primitive
is true
(default) the primitive unit cell is returned, else the conventional unit cell is returned.
DFTK.symmetries_preserving_kgrid
— MethodFilter out the symmetry operations that respect the symmetries of the discrete BZ grid
DFTK.symmetrize
— MethodSymmetrize a RealFourierArray
by applying all the model symmetries (by default) and forming the average.
DFTK.symmetry_operations
— FunctionReturn the $k$-point symmetry operations associated to a lattice and atoms.
DFTK.total_local_potential
— MethodGet the total local potential of the given Hamiltonian, in real space.
DFTK.unit_to_au
— Methodunit_to_ao(symbol)
Get the factor converting from the unit symbol
to atomic units. E.g. unit_to_au(:eV)
returns the conversion factor from electron volts to Hartree.
DFTK.weighted_ksum
— MethodSum an array over kpoints, taking weights into account
DFTK.ylm_real
— MethodReturns the (l,m) real spherical harmonic Ylm(r). Consistent with https://en.wikipedia.org/wiki/Tableofsphericalharmonics#Realsphericalharmonics
DFTK.@timing
— MacroShortened version of the @timeit
macro from TimerOutputs
, which writes to the DFTK timer.
DFTK.@timing_seq
— MacroSimilar to @timing
, but disabled in parallel runs. Should be used to time threaded regions, since TimerOutputs is not thread-safe and breaks otherwise.
DFTK.Smearing.entropy
— MethodEntropy. Note that this is a function of the energy x
, not of occupation(x)
. This function satisfies s' = x f' (see https://www.vasp.at/vasp-workshop/k-points.pdf p. 12 and https://arxiv.org/pdf/1805.07144.pdf p. 18.
DFTK.Smearing.occupation
— MethodOccupation at x
, where in practice x = (ε - εF) / T.
DFTK.Smearing.occupation_derivative
— MethodDerivative of the occupation function, approximation to minus the delta function.
DFTK.Smearing.occupation_divided_difference
— Method(f(x) - f(y))/(x - y), computed stably in the case where x and y are close