DFWannier.BVectorShellsType
BVectorShells

Shells of bvectors.

The bvectors are sorted by norm such that equal-norm bvectors are grouped into one shell.

Fields

  • recip_lattice: 3 * 3, each column is a reciprocal lattice vector
  • kpoints: 3 * n_kpts, in fractional coordinates
  • bvectors: vectors of 3 * n_bvecs_per_shell, in cartesian coordinates
  • weights: vector of float, weights of each shell
  • multiplicities: number of bvectors in each shell
  • n_bvecs: total number of bvectors
DFWannier.BVectorShellsMethod
BVectorShells(recip_lattice, kpoints, bvectors, weights)

Constructor of BVectorShells.

Only essential arguments are required, remaing fields of BVectorShells are initialized accordingly. This should be used instead of directly constructing BVectorShells.

Arguments

  • recip_lattice: 3 * 3, each column is a reciprocal lattice vector
  • kpoints: 3 * n_kpts, in fractional coordinates
  • bvectors: vectors of 3 * n_bvecs_per_shell, in cartesian coordinates
  • weights: vector of float, weights of each shell
DFWannier.BVectorsType
BVectors

The bvectors for each kpoint.

Fields

  • recip_lattice: 3 * 3, each column is a reciprocal lattice vector
  • kpoints: 3 * n_kpts, in fractional coordinates
  • bvectors: 3 * n_bvecs, in cartesian coordinates
  • weights: n_bvecs, weights of each bvector
  • kpb_k: k+b vectors at kpoint k, k -> k + b (index of periodically equivalent kpoint inside recip_lattice)
  • kpb_b: 3 * n_bvecs * n_kpts, displacements between k + b and its periodic image inside recip_lattice, such that k+b = kpoints[:, kpb_k[ib, ik]] + kpb_b[:, ib, ik] (in fractional)
  • n_kpts: number of kpoints
  • n_bvecs: total number of bvectors
Note

In principle, we don't need to sort the bvectors for each kpoint, so that the bvectors have the same order as each kpoint. However, since Wannier90 sort the bvectors, and the mmn file is written in that order, so we also sort in the same order as Wannier90.

DFWannier.ColinMatrixType
ColinMatrix{T, M <: AbstractMatrix{T}} <: AbstractMagneticMatrix{T}

Defines a Hamiltonian Matrix with [up zeros; zeros down] structure. It is internally only storing the up and down block.

DFWannier.Exchange2ndOrderType
Exchange2ndOrder{T <: AbstractFloat}

This holds the exhanges between different orbitals and calculated sites. Projections and atom datablocks are to be found in the corresponding wannier input file. It turns out the ordering is first projections, then atom order in the atoms datablock.

DFWannier.Exchange4thOrderType
Exchange4thOrder{T <: AbstractFloat}

This holds the exhanges between different orbitals and calculated sites. Projections and atom datablocks are to be found in the corresponding wannier input file. It turns out the ordering is first projections, then atom order in the atoms datablock.

DFWannier.HamiltonianKGridMethod
HamiltonianKGrid(hami::TBHamiltonian{T}, nk, H_function_k::Function = x -> nothing) where T
HamiltonianKGrid(hami::TBHamiltonian{T}, k_grid, H_function_k::Function = x -> nothing) where T

Takes a k grid, calculates Hk for each of them and diagonalizes. Only the eigenvectors and eigenvalues of Hk are stored, the H_function_k function is called on the intermediate Hk.

DFWannier.MagneticVectorType

Vector following the same convention as the in AbstractMagneticMatrix, i.e. first half of the indices contain the up part, second the down part

DFWannier.NonColinMatrixType
NonColinMatrix{T, M <: AbstractMatrix{T}} <: AbstractMagneticMatrix{T}

Defines a Hamiltonian Matrix with [up updown;downup down] structure. Since atomic projections w.r.t spins are defined rather awkwardly in Wannier90 for exchange calculations, i.e. storing the up-down parts of an atom sequentially, a NonColinMatrix reshuffles the entries of a matrix such that it follows the above structure.

DFWannier.TBBlockType
TBBlock

Building block for TBOperator. It holds the matrix elements of the operator between central and a shifted unit cell. Upon construction, the wigner-seitz shifts are taken into account to create the correct matrix elements between the Wannierfunctions, stored in tb_block. The block field is basically tb_block but with each element divided by the amount of Wigner-Seitz degeneracies and shifts which speeds up later k-point interpolation.

DFWannier.TBOperatorType
TBOperator

Alias for a Vector of TBBlocks. Indexing with NTuple{3,Int} or Vec3 is supported which allows for easily retrieving the TBBlock that corresponds to the shifted unit cell. Aliases: TBHamiltonian, TBSpin

Base.convertMethod

Reshuffles standard Wannier90 up-down indices to the ones for the structure of a NonColinMatrix.

Base.convertMethod

Reshuffles standard Wannier90 up-down indices to the ones for the structure of a MagneticVector.

DFWannier.DHvecvalsMethod
DHvecvals(hami::TBHamiltonian{T, Matrix{T}}, k_grid::Vector{Vec3{T}}, atoms::Atom{T}) where T <: AbstractFloat

Calculates $D(k) = [H(k), J]$, $P(k)$ and $L(k)$ where $H(k) = P(k) L(k) P^{-1}(k)$. hami should be the full Hamiltonian containing both spin-diagonal and off-diagonal blocks.

DFWannier.HkMethod
Hk(hamiltonian::TBHamiltonian, kpoint::Vec3)
Hk!(hk::AbstractMatrix, hamiltonian::TBHamiltonian, kpoint::Vec3)

Constructs the reciprocal Hamiltonian at a given k-point.

DFWannier.S_RMethod
S_R(chk, Sx, Sy, Sz)

Takes the DFT Sx, Sy, Sz spin matrices and constructs the TBSpin from them. Using the Wannier90 checkpoint information in chk.

DFWannier._bvec_to_kbMethod

Find equivalent kpoint and displacement vector of bvectors bvecs at kpoint k.

all inputs in fractional coordinates.

DFWannier._sort_kbMethod

Sort bvectors specified by equivalent kpoint indices k and cell displacements b.

Sorting order:

  1. length of bvectors: nearest k+b goes first, this is achieved by comparing the norm bvecs_norm.
  2. supercell index: the supercell are already sorted by sort_supercell, which generates our input translations.
  3. index of kpoint: the smaller index goes first, dictated by the input kpoints.

bvecs_norm: length of bvectors, cartesian norm. k: index in kpoints for equivalent kpoint of bvectors b: cell displacements of bvectors, fractional coordinates translations: of supercell, fractional coordinates

DFWannier.are_parallelMethod
are_parallel(A, B; atol=1e-6)

Check if the columns of matrix A and columns of matrix B are parallel.

Arguments

  • A: matrix
  • B: matrix

Keyword Arguments

  • atol: tolerance to check parallelism
DFWannier.calc_angmomMethod

Calculates the angular momentum between two wavefunctions and around the center.

DFWannier.calc_dipMethod

Calculates the dipole term between two wavefunctions. Make sure the wavefunctions are normalized!

DFWannier.calc_exchangesMethod
calc_exchanges(hamiltonian::TBHamiltonian, atoms::Vector{<:Atom}, fermi, exchange_type; kwargs...)

Calculates the magnetic exchange parameters between the atoms. exchange_type can be Exchange2ndOrder or Exchange4thOrder. The kwargs control various numerical parameters for the calculation:

  • nk = (10,10,10): the amount of k-points to be used for the uniform interpolation grid.
  • R = (0,0,0): the unit cell index to which the exchange parameters are calculated.
  • ωh = -30.0: the lower bound of the energy integration
  • ωv = 0.15: the height of the contour in complex space to integrate the Green's functions
  • n_ωh = 3000: number of integration points along the horizontal contour direction
  • n_ωv = 500: number of integration points along the vertical contour direction
  • site_diagonal = false: if true the hamiltonians and Δ will diagonalized on-site and the

returned exchange matrices hold the exchanges between well-defined orbitals. If this is not done, the exchange matrices entries don't mean anything on themselves and a trace should be performed to find the exchange between the spins on sites i and j.

DFWannier.check_b1Method
check_b1(shells::BVectorShells; atol=1e-6)

Check completeness (B1 condition) of BVectorShells.

Arguments

  • shells: BVectorShells containing bvectors in each shell

Keyword Arguments

  • atol: tolerance, equivalent to Wannier90 input parameter kmesh_tol
DFWannier.compute_weightsMethod
compute_weights(bvectors::Vector{Matrix{T}}; atol=1e-6)

Try to guess bvector weights from MV1997 Eq. (B1).

The input bvectors are overcomplete vectors found during shell search, i.e. from search_shells. This function tries to find the minimum number of bvector shells that satisfy the B1 condition, and return the new BVectorShells and weights.

Arguments

  • bvectors: vector of bvectors in each shell

Keyword Arguments

  • atol: tolerance to satisfy B1 condition, equivalent to Wannier90 input parameter kmesh_tol
DFWannier.compute_weightsMethod
compute_weights(shells::BVectorShells; atol=1e-6)

Try to guess bvector weights from MV1997 Eq. (B1).

Arguments

  • shells: BVectorShells containing bvectors in each shell

Keyword Arguments

  • atol: tolerance to satisfy B1 condition, equivalent to Wannier90 input parameter kmesh_tol
DFWannier.delete_parallelMethod
delete_parallel(shells::BVectorShells)

Remove shells having parallel bvectors.

Arguments

  • shells: BVectorShells containing bvectors in each shell
DFWannier.delete_parallelMethod
delete_parallel(bvectors::Vector{Matrix{T}})

Remove shells having parallel bvectors.

Arguments

  • bvectors: vector of bvectors in each shell
DFWannier.flatten_shellsMethod
flatten_shells(shells::BVectorShells)

Flatten shell vectors into a matrix.

Return a tuple of (bvecs, bvecs_weight), where

  • bvecs: 3 * n_bvecs
  • bvecs_weight: n_bvecs
DFWannier.fourier_q_to_RMethod
fourier_q_to_R(f::Function, q_vectors, R_vectors)

Performs a fourier transform from the ab-initio kpoints to the wigner seitz unit cells. The function will be executed inside the fourier transform loop, being called like f(iR, ik, phase)

DFWannier.fourier_transformMethod

Fourier transforms the tight binding hamiltonian and calls the R_function with the current index and the phase.

DFWannier.generate_TBBlocksMethod
generate_TBBlocks(chk::NamedTuple, O_R::Vector)

Generates the Vector of TBBlocks from the Wannier90 checkpoint info in chk, and the real space operator O_R. This preapplies the Wigner Seitz shifts and degeneracies to speed up k-point interpolation.

DFWannier.get_bvectorsMethod
get_bvectors(kpoints, recip_lattice; kmesh_tol=1e-6)

Generate and sort bvectors for all the kpoints.

Arguments

  • kpoints: 3 * n_kpts, kpoints in fractional coordinates
  • recip_lattice: 3 * 3, columns are reciprocal lattice vectors

Keyword Arguments

  • kmesh_tol: equivalent to Wannier90 input parameter kmesh_tol
DFWannier.get_bvectors_nearestMethod
get_bvectors_nearest(kpoints, recip_lattice; kmesh_tol=1e-6)

Generate and sort bvectors for all the kpoints.

Arguments

  • kpoints: 3 * n_kpts, kpoints in fractional coordinates
  • recip_lattice: 3 * 3, columns are reciprocal lattice vectors

Keyword Arguments

  • kmesh_tol: equivalent to Wannier90 input parameter kmesh_tol
DFWannier.index_bvectorMethod
index_bvector(kpb_k, kpb_b, k1, k2, b)

Given bvector b connecting kpoints k1 and k2, return the index of the bvector ib.

This is a reverse search of bvector index if you only know the two kpoints k1 and k2, and the connecting displacement vector b.

Arguments

  • kpb_k: n_bvecs * n_kpts, k+b kpoints at k1
  • kpb_b: 3 * n_bvecs * n_kpts, displacement vector for k+b bvectors at k1
  • k1: integer, index of kpoint k1
  • k2: integer, index of kpoint k2
  • b: vector of 3 integer, displacement vector from k1 to k2
DFWannier.make_supercellFunction
make_supercell(kpoints, replica=5)

Make a supercell of kpoints by translating it along 3 directions.

Arguments

  • replica: integer, number of repetitions along ±x, ±y, ±z directions
DFWannier.make_supercellMethod
make_supercell(kpoints, replica)

Make a supercell of kpoints by translating it along 3 directions.

On output there are (2*replica + 1)^3 cells, in fractional coordinates.

Arguments

  • kpoints: 3 * n_kpts, in fractional coordinates
  • replica: 3, number of repetitions along ±x, ±y, ±z directions
DFWannier.r_RMethod
r_R(chk, kbonds)

Constructs the r TBOperator from the Wannier90 checkpoint info chk and the kbond information that can be read with read_nnkp.

DFWannier.read_chkMethod
read_chk(chk_file)
read_chk(job::Job)

Reads a Wannier90 .chk file returning a NamedTuple containing all the information.

DFWannier.read_colin_hamiMethod
read_colin_hami(up_chk, down_chk, up_eigvals::AbstractString, down_eigvals::AbstractString)

Returns the colinear TBHamiltonian representing the up-down blocks of the Wannier Tight Binding Hamiltonian.

DFWannier.read_eigMethod
read_eig(eig_file)

Reads the DFT eigenvalues from a .eig file.

DFWannier.read_hamiltonianMethod
read_hamiltonian(job::Job)

Goes through the job and will attempt to read the hamiltonian files. If it finds a colinear calculation in the job it will read the up and down hamiltonians, if the job was either nonmagnetic or noncolinear it will read only one hamiltonian file (there should be only one).

DFWannier.read_hamiltonianMethod
read_hamiltonian(chk::NamedTuple, eigvals::Matrix)

Uses the Wannier90 chkpoint info in chk and DFT eigenvals read with [read_eig] to construct the TBHamiltonian.

DFWannier.read_nnkpMethod
read_nnkp(nnkp_file)

Reads a Wannier90 .nnkp file and returns (recip_cell, kpoints, kbonds).

DFWannier.read_rMethod
read_r(chk_file::AbstractString, nnkp_file::AbstractString)
read_r(job::Job)

Constructs the r [TBOperator] from the Wannier90 .chk and .nnkp files. This requires that the k_neighbor_weights is written into the .chk file and might need a patched Wannier90 version.

DFWannier.read_spinMethod
read_spin(chk_file, spn_file)
read_spin(job::Job)

Reads the .spn and .chk files to generate a TBSpin tight-binding spin operator.

DFWannier.read_spnMethod
read_spn(filename)

Reads a .spn file and returns the DFT Sx, Sy, Sz. They are a Vectors with nk Matrices of size (nb, nb), where nk is the number of k-points and nb the number of bands.

DFWannier.read_wannier_blocksMethod
read_wannier_blocks(f)

Reads a Wannier90 file such as .nnkp and separates each begin end block into an entry in the ouput Dict.

DFWannier.search_shellsMethod
search_shells(kpoints, recip_lattice; atol=1e-6, max_shells=36)

Search bvector shells satisfing B1 condition.

Arguments

  • kpoints: fractional coordinates
  • recip_lattice: each column is a reciprocal lattice vector

Keyword Arguments

  • atol: tolerance to select a shell (points having equal distances), equivalent to Wannier90 input parameter kmesh_tol.
  • max_shells: max number of nearest-neighbor shells, equivalent to Wannier90 input parameter search_shells.
DFWannier.sort_bvectorsMethod
sort_bvectors(shells::BVectorShells; atol=1e-6)

Sort bvectors in shells at each kpoints, to be consistent with Wannier90.

Wannier90 use different order of bvectors at each kpoint, in principle, this is not needed. However, the mmn file is written in such order, so we need to sort bvectors and calculate weights, since nnkp file has no section of weights.

Arguments

  • shells: BVectorShells

Keyword Arguments

  • atol: equivalent to Wannier90 input parameter kmesh_tol
DFWannier.sort_supercellMethod
sort_supercell(translations, recip_lattice; atol=1e-8)

Sort supercell to fix the order of bvectors.

Both input and output translations are in fractional coordinates.

Arguments

  • translations: 3 * n_supercell matrix, in fractional coordinates
  • recip_lattice: each column is a reciprocal lattice vector

Keyword Arguments

  • atol: tolerance to compare bvectors, this is the same as what is hardcoded in Wannier90
Note

This is used to reproduce Wannier90 bvector order.

DFWannier.wannierbandsMethod
wannierbands(hamiltonian::TBHamiltonian, kpoints::Vector{Vec3})
wannierbands(hamiltonian::TBHamiltonian, bands::Vector{DFControl.AbstractBand}

Constructs the whole bandstructure for a given set of k-points and TBHamiltonian.

DFWannier.write_xsfMethod
write_xsf(filename::String, wfc::WannierFunction, str::Structure; value_func=x->norm(x))

Writes a WannierFunction and Structure to an xsf file that is readable by XCrysden or VESTA. The values that are written can be controlled by value_func that gets used on each entry of wfc.values and should output a single Number.

DFWannier.σxMethod

Generates a Pauli σx matrix with the dimension that is passed through n.

DFWannier.σyMethod

Generates a Pauli σy matrix with the dimension that is passed through n.

DFWannier.σzMethod

Generates a Pauli σz matrix with the dimension that is passed through n.