Crystal data

General

Basic data grids

Electrum.DataGridType
Electrum.DataGrid{D,B<:LatticeBasis,S<:AbstractVector{<:Real},T} <: AbstractArray{T,D}

Stores a grid of values of type T defined in a D-dimensional crystal lattice basis of type B with a shift parameter of type S.

By convention, indexing of DataGrid is zero-based. This convention is used so that the first entry corresponds to data at the origin can be indexed with zeros. However, getindex() is implemented such that the dataset may be indexed by any integer, with modulo math used to convert to an index within the grid. This is done with Electrum.reinterpret_index().

Linear indexing is one-based like that of the underlying array.

Type aliases

For convenience, the aliases RealDataGrid and ReciprocalDataGrid are provided and are defined below:

const RealDataGrid{D} = DataGrid{RealBasis{D,Float64},SVector{D,Float64},D}
const ReciprocalDataGrid{D} = DataGrid{ReciprocalBasis{D,Float64},KPoint{D},D}

Note that ReciprocalDataGrid uses a KPoint{D} to represent the shift, as opposed to an SVector{D}. This allows for the concurrent storage of a weight along with the shift, which may be relevant for wavefunctions which exploit the symmetry of the k-point mesh.

Electrum.RealDataGridType
Electrum.DataGrid{D,B<:LatticeBasis,S<:AbstractVector{<:Real},T} <: AbstractArray{T,D}

Stores a grid of values of type T defined in a D-dimensional crystal lattice basis of type B with a shift parameter of type S.

By convention, indexing of DataGrid is zero-based. This convention is used so that the first entry corresponds to data at the origin can be indexed with zeros. However, getindex() is implemented such that the dataset may be indexed by any integer, with modulo math used to convert to an index within the grid. This is done with Electrum.reinterpret_index().

Linear indexing is one-based like that of the underlying array.

Type aliases

For convenience, the aliases RealDataGrid and ReciprocalDataGrid are provided and are defined below:

const RealDataGrid{D} = DataGrid{RealBasis{D,Float64},SVector{D,Float64},D}
const ReciprocalDataGrid{D} = DataGrid{ReciprocalBasis{D,Float64},KPoint{D},D}

Note that ReciprocalDataGrid uses a KPoint{D} to represent the shift, as opposed to an SVector{D}. This allows for the concurrent storage of a weight along with the shift, which may be relevant for wavefunctions which exploit the symmetry of the k-point mesh.

Electrum.ReciprocalDataGridType
Electrum.DataGrid{D,B<:LatticeBasis,S<:AbstractVector{<:Real},T} <: AbstractArray{T,D}

Stores a grid of values of type T defined in a D-dimensional crystal lattice basis of type B with a shift parameter of type S.

By convention, indexing of DataGrid is zero-based. This convention is used so that the first entry corresponds to data at the origin can be indexed with zeros. However, getindex() is implemented such that the dataset may be indexed by any integer, with modulo math used to convert to an index within the grid. This is done with Electrum.reinterpret_index().

Linear indexing is one-based like that of the underlying array.

Type aliases

For convenience, the aliases RealDataGrid and ReciprocalDataGrid are provided and are defined below:

const RealDataGrid{D} = DataGrid{RealBasis{D,Float64},SVector{D,Float64},D}
const ReciprocalDataGrid{D} = DataGrid{ReciprocalBasis{D,Float64},KPoint{D},D}

Note that ReciprocalDataGrid uses a KPoint{D} to represent the shift, as opposed to an SVector{D}. This allows for the concurrent storage of a weight along with the shift, which may be relevant for wavefunctions which exploit the symmetry of the k-point mesh.

AbstractFFTs.fftFunction
fft(g::RealDataGrid) -> ReciprocalDataGrid
fft(g::ReciprocalDataGrid) -> RealDataGrid

Performs a fast Fourier transform on the data in a DataGrid.

AbstractFFTs.ifftFunction
ifft(g::RealDataGrid) -> ReciprocalDataGrid
ifft(g::ReciprocalDataGrid) -> RealDataGrid

Performs an inverse fast Fourier transform on the data in a DataGrid.

The inverse FFT is normalized so that ifft(fft(g)) ≈ g (to within floating point error).

AbstractFFTs.fftfreqMethod
fftfreq(g::DataGrid{D}) -> Array{SVector{D,Float64},D}

Returns the Fourier transform frequency bins for an DataGrid.

For real space data, the frequency bins will be angular wavenumbers, matching the 2π factors that are introduced when transforming between a RealBasis and a ReciprocalBasis. The convention used by FFTW.fftfreq() is also used: frequency bins at or above the Nyquist frequency will be negative.

For reciprocal space data, the frequencies are binned with the assumption that the lattice vectors are given in angular wavenumbers, and they represent real space coordinates. The Nyquist frequency convention is not used, so all elements will have positive indices.

Electrum.volumeMethod
volume(g::DataGrid) -> eltype(basis(g))

Gets the crystal volume associated with a DataGrid. Units are assumed to be bohr³ for DataGrid types with a RealBasis, and rad³*bohr⁻³ for types with a ReciprocalBasis.

Electrum.voxelsizeFunction
voxelsize(g::DataGrid) -> eltype(basis(g))

Gets the size of a single real space voxel associated with a data grid. Units are assumed to be bohr³.

Electrum.integrateFunction
integrate([f::Function = identity], g::RealDataGrid{D,T}) -> T

Applies the function f elementwise to a datagrid, then integrates the grid across all voxels, accounting for the voxel volumes.

Electrum.partial_derivativeFunction
partial_derivative(g::RealDataGrid, dim::Integer)

Calculates the partial derivative of g along coordinates with respect to the basis vector of dimension dim.

Electrum.cell_gradientFunction
cell_gradient(g::RealDataGrid{D}) -> RealDataGrid{D,<:SVector{D}}

Calculates the gradient of a RealDataGrid, which consists of vectors containing the components of the partial derivatives along coordinates with respect to the basis vectors of the unit cell. The gradient may need to be transformed to obtain a more standard result: for this, use the gradient() function.

Electrum.gradientFunction
gradient(g::RealDataGrid{D}) -> RealDataGrid{D,<:SVector{D}}

Calculates the gradient of a RealDataGrid, which consists of vectors containing the components of the partial derivatives along the orthonormal real space basis. The units of the output are those of the input multiplied by inverse bohr.

Electrum.directional_derivativeFunction
directional_derivative(g::RealDataGrid{D}, v::AbstractVector) -> RealDataGrid{D}

Calculates the directional derivative of a RealDataGrid along a vector v in Cartesian coordinates. Note that v does not need to be normalized, nor is it automatically normalized; an unnormalized input will scale the results accordingly.

Electrum.remove_shiftFunction
remove_shift(g::RealDataGrid) -> RealDataGrid

Shifts the data in a RealDataGrid so that the shift vector of the new RealDataGrid is zero. Currently, only shifts for trivial cases (where the data shifting can be accomplished by a simple circular shift of the array) are implemented; future implementations will need to use Fourier transforms to accomplish this task.

FFT iterators

Electrum.FFTLengthType
FFTLength <: AbstractVector{Int}

The one-dimensional counterpart to FFTBins, supporting only a single dimension. Its elements are plain Int types rather than the CartesianIndex of FFTBins.

In essence, it serves as a counterpart to Base.OneTo for FFT bins.

Examples

julia> Electrum.FFTLength(4)
4-element FFTLength:
 0
 1
 -2
 -1
Electrum.FFTBinsType
FFTBins{D} <: AbstractArray{CartesianIndex{D},D}

An iterable object defining a range of integer FFT bins. This can be used to convert a Cartesian array index to an FFT bin index, or vice versa.

The outputs use the convention where frequencies at or above the Nyquist frequency for that dimension are negative, matching the output of FFTW.fftfreq.

Examples

julia> FFTBins(4)
4-element FFTIndices{1}:
 CartesianIndex(0,)
 CartesianIndex(1,)
 CartesianIndex(-2,)
 CartesianIndex(-1,)

julia> FFTBins(3, 3)
3×3 FFTIndices{2}:
 CartesianIndex(0, 0)   CartesianIndex(0, 1)   CartesianIndex(0, -1)
 CartesianIndex(1, 0)   CartesianIndex(1, 1)   CartesianIndex(1, -1)
 CartesianIndex(-1, 0)  CartesianIndex(-1, 1)  CartesianIndex(-1, -1)

k-points

Electrum.KPointType
KPoint{D} <: DenseVector{Float64}

Stores a k-point with an associated weight that corresponds to the number of symmetry-equivalent k-points, stored as an integer.

Electrum.KPointMeshType
KPointMesh{D} <: AbstractVector{KPoint{D}}

Contains a list of k-points associated with a matrix describing the mesh that was used to generate the points, and its shift off the Γ point (origin). If the mesh used to generate the points is unknown, it will be set to the zero matrix of dimension D.

A KPointMesh can be indexed as if it were an ordinary Vector{KPoint{D}}.

Electrum.nkptFunction
nkpt(k::KPointMesh) -> Int

Counts the number of k-points (equivalent to length(k)). This function can be defined for custom types that contain a KPointMesh.

nkpt(x) -> Int

Counts the number of explicitly enumerated k-points associated with an object containing a KPointMesh.

By default, this function returns length(KPointMesh(x)), so defining KPointMesh(x::T) for a type T containing a KPointMesh will automatically define nkpt(x::T).

Electrum.weightFunction
weight(k::KPoint) -> Int

Returns the weight associated with a k-point.

Spin states

Electrum.MultiplicityType
Multiplicity{M} <: AbstractUnitRange{Rational{Int}}

Represents a valid range of spin states corresponding with multiplicity M. These index as if they are UnitRange{Rational{Int}} objects of the form -(M - 1)//2:(M - 1)//2, and can be converted to UnitRange objects with the UnitRange constructor.

These types are singleton types with a numeric Int parameter, allowing for its use as a type parameter in other types which require information on the spin multiplicity or allowed spin states.

Examples

julia> UnitRange(SpinRange{3}())
1:1

julia> SpinRange{4}() == -3//2:3//2
true
Electrum.SpinBivectorType
SpinBivector{D,T}

A skew-symmetric matrix representing a spin bivector in D dimensions. No automatic normalization has been implemented for this type.

Why the spin axis is a bivector

The identification of a spin direction with a vector is only possible in 3D. This is because the axis of rotation is defined by the Hodge dual of the rotation plane. The Hodge dual is a linear map relating k-dimensional objects in D-dimensional space to D-k-dimensional objects in the same space.

A bivector is a 2-dimensional object, and in three dimensions, its Hodge dual is a vector. This is not the case in any other number of dimensions. To allow for generality in describing spin, this representation is used instead.

Energies and occupancies

Electrum.EnergyOccupancyType
EnergyOccupancy{T<:Real}

A data structure consisting of a pair of an energy value and an occupancy number, both of type T. Energies are assumed to be in Hartree. Occupancy values are not constrained, but will generally range from 0 to 2 for the results of restricted calculations (no separate treatment of spins), or from 0 to 1 for unrestricted calculations (wavefunctions with separate spins).

Electrum.EnergiesOccupanciesType
EnergiesOccupancies{T,N} <: AbstractArray{EnergyOccupancy{T},N}

Type alias for Array{EnergyOccupancy{T},N}. Data structures S which contain EnergyOccupancy values in a collection should define the constructor EnergiesOccupancies(::S).

Electrum.energyFunction
energy(eo::EnergyOccupancy{T}) -> T

Returns the energy value in an EnergyOccupancy.

Electrum.occupancyFunction
occupancy(eo::EnergyOccupancy{T}) -> T

Returns the occupancy value in an EnergyOccupancy.

Electrum.energiesFunction
energies(x) -> Array{<:Real}

Returns the energy data associated with a collection of EnergyOccupancy{T} objects. By default, this falls back to energy.(EnergyOccupancy(x)).

energies(d::AbstractDensityOfStates; usefermi=false) -> Vector{Float64}

Gets the range of energies in the dataset. If usefermi is set to true, the energies returned will be adjusted such that the Fermi energy is set to zero.

There are no guarantees on the unit of energy used!

Electrum.occupanciesFunction
occupancies(x) -> Array{<:Real}

Returns the occupancy data associated with a collection of EnergyOccupancy{T} objects. By default, this falls back to occupancy.(EnergyOccupancy(x)).

Electrum.min_energyFunction
min_energy(x) -> Real

Returns the minimum energy in a collection of EnergyOccupancy data or an object containing such data.

Electrum.max_energyFunction
max_energy(x) -> Real

Returns the maximum energy in a collection of EnergyOccupancy data or an object containing such data.

Electrum.min_occupancyFunction
min_occupancy(x) -> Real

Returns the minimum occupancy in a collection of EnergyOccupancy data or an object containing such data. For most raw calculation data, this should return zero if unoccupied states were considered.

Electrum.max_occupancyFunction
max_occupancy(x) -> Real

Returns the maximum occupancy in a collection of EnergyOccupancy data or an object containing such data. For a restricted calculation (no explicit treatment of spin), this is usually around 2, and for an unrestricted calculation (explicit spin treatment) this is usually around 1.

In many cases, you may want to determine the maximum possible occupancy value, not the maximum in the dataset, in which case, you should use round(Int, max_occupancy(a), RoundUp).

Electrum.fermiFunction
fermi(x) -> Real

Estimates the Fermi energy using the provided energy and occupancy data by interpolating the data between the occupancies nearest half of the maximum occupancy.

fermi(d::AbstractDensityOfStates) -> Float64

Gets the Fermi energy from DOS data. There are no guarantees on the unit of energy used!

Wavefunctions

Electrum.PlanewaveWavefunctionType
PlanewaveWavefunction{D,T} <: AbstractArray{T,D}

Stores the components of a wavefunction constructed from a planewave basis. Usually, the coefficient data type T will be a ComplexF32, as in DFT calculations, double precision convergence of the density will correspond to single-precision converegnce of the wavefunction.

Internally, coefficients are stored in an Array{4,T}. Indexing is then manually implemented, with a D-dimensional CartesianIndex used for accessing each coefficient associated with a G-vector. PlanewaveWavefunction instances are mutable, with getindex() and setindex!() defined for them, but they are not resizable, and the backing array should not be resized.

Electrum.PlanewaveIndexType
Electrum.PlanewaveIndex{D}

A special indexing type used to index the components of wavefunctions in a planewave basis.

In many computational chemistry packages, the standard indexing of wavefunction components occurs in the following canonical order: spins, k-points, bands, then the h, k, and l indices of the G-vectors associated with the coefficients. However, in many cases, users will want to select a spin, k-point, or band before selecting a G-vector index.

To keep the syntax intuitive while maintaining performance (mostly by ensuring that all of the components of a wavefunction are stored compactly) this index type ensures that the iteration occurs in a natural order for users and in an efficient order for Julia.

This also ensures that G-vectors with negative indices are handled correctly: while the canonical G-vectors are those within some defined ranges of indices, out of bounds indices are reinterpreted automatically using modulo arithmetic.

Electrum.nspinFunction
nspin(wf::PlanewaveWavefunction) -> Int

Returns the number of spins associated with a PlanewaveWavefunction.

Electrum.nbandMethod
nband(wf::PlanewaveWavefunction) -> Int

Returns the number of bands (occupied and unoccupied) associated with a PlanewaveWavefunction.

Electrum.nonzero_g_indicesFunction
nonzero_g_indices(wf::PlanewaveWavefunction{D}) -> Vector{CartesianIndex{D}}

Returns a vector of CartesianIndex objects corresponding to planewave G-vector indices that are not all zero at each band and k-point.

To return the G-vectors as objects which subtype AbstractVector, nonzero_g_vectors(wf) may be used instead, which is equivalent to calling SVector.(Tuple.(nonzero_g_indices(wf))).

Electrum.nonzero_g_vectorsFunction
nonzero_g_vectors(wf::PlanewaveWavefunction{D}) -> Vector{SVector{D,Int}}

Returns a vector of SVector{D,Int} objects corresponding to planewave G-vectors that are not all zero at each band and k-point.

This is equivalent to calling SVector.(Tuple.(nonzero_g_indices(wf))), and can be used whenever AbstractVector inputs are needed instead of a CartesianIndex.

Band structures

Electrum.BandAtKPointType
BandAtKPoint

Stores information about a band's energy and its occupancy at a specific k-point.

Electrum.BandStructureType
BandStructure{D}

Stores information about an electronic band structure, including the list of k-points used to generate the data (as an AbstractVector{KPoint{D}}) and the band information at every k-point (as a Vector{BandAtKPoint}).

Electrum.nbandMethod
nband(b::BandStructure) -> Int

Returns the number of bands (occupied and unoccupied) associated with a BandStructure.

Electrum.nbandMethod
nband(b::BandAtKPoint) -> Int

Returns the number of bands associated with a k-point.

Electrum.FatBandsType
FatBands{D}

Stores information relevant to plotting fatbands.

  • FatBands.bands: matrix of energies at each [kpt, band].
  • FatBands.projband: array of lm-decomposed band structure. [orbital, ion, band, kpt].
  • FatBands.cband: array of complex-valued contributions to band structure.

Density of states

Electrum.fermiMethod
fermi(d::AbstractDensityOfStates) -> Float64

Gets the Fermi energy from DOS data. There are no guarantees on the unit of energy used!

Electrum.smearFunction
smear(dos::DensityOfStates, sigma::Real)

Smears the DOS function by convoluting it with a normalized Gaussian function with width sigma.

Electrum.energiesMethod
energies(d::AbstractDensityOfStates; usefermi=false) -> Vector{Float64}

Gets the range of energies in the dataset. If usefermi is set to true, the energies returned will be adjusted such that the Fermi energy is set to zero.

There are no guarantees on the unit of energy used!

Electrum.nelectronsFunction
nelectrons(d::DensityOfStates)

Gets the approximate number of electrons that are needed to reach the Fermi level.

Atomic data

Electrum.SphericalHarmonicType
SphericalHarmonic{Lmax}

Real spherical harmonic components up to Lmax. This can be used to describe atomic orbitals or projections of data onto atomic sites.

LinearAlgebra.dotMethod
dot(sh1::SphericalHarmonics, sh2::SphericalHarmonics) -> Float64

Calculates the dot product between the components of two spherical harmonics, which can be used to measure the degree of similarity between them. Note that this does not account for differences in rotation between the spherical harmonics.