Initializing a System

The first step in the wannierberri calculation is initialising the System. This is done by means of child classes System described below. They all have an important common method set_symmetry(). The system may come either from Wanier functions constructed by Wannier90, or from ref:tight binding <sec-tb-models> models.

class wannierberri.system.System[source]

Bases: object

The base class for describing a system. Does not have its own constructor, please use the child classes, e.g System_w90 or System_tb

Parameters:
  • seedname (str) – the seedname used in Wannier90. Default: wannier90

  • berry (bool) – set True to enable evaluation of external term in Berry connection or Berry curvature and their derivatives. Default: False

  • spin (bool) – set True if quantities derived from spin will be used. Default:False

  • morb (bool) – set True to enable calculation of external terms in orbital moment and its derivatives. Requires the .uHu file. Default: False

  • periodic ([bool,bool,bool]) – set True for periodic directions and False for confined (e.g. slab direction for 2D systems). If less then 3 values provided, the rest are treated as False . Default : (True, True, True)

  • SHCryoo (bool) – set True if quantities derived from Ryoo’s spin-current elements will be used. (RPS 2019) Default: False

  • SHCqiao (bool) – set True if quantities derived from Qiao’s approximated spin-current elements will be used. (QZYZ 2018). Default: False

  • use_ws (bool) – minimal distance replica selection method Minimal-distance replica selection method. equivalent of use_ws_distance in Wannier90. Default: True

  • mp_grid ([nk1,nk2,nk3]) – size of Monkhorst-Pack frid used in ab initio calculation. Needed when use_ws=True, and only if it cannot be read from input file, i.e. like System_tb, System_PythTB, System_TBmodels ,:class:.System_fplo, but only if the data originate from ab initio data, not from toy models. In contrast, for System_w90 and System_ase it is not needed, but can be provided and will override the original value (if you know what and why you are doing) Default: None

  • frozen_max (float) – position of the upper edge of the frozen window. Used in the evaluation of orbital moment. But not necessary. Default: -inf

  • _getFF (bool) – generate the FF_R matrix based on the uIu file. May be used for only testing so far. Default : False

  • use_wcc_phase (bool) – using wannier centers in Fourier transform. Correspoinding to Convention I (True), II (False) in Ref.”Tight-binding formalism in the context of the PythTB package”. Default: False

  • wannier_centers_cart (array-like(num_wann,3)) – use the given wannier_centers (cartesian) instead of those determined automatically. Incompatible with wannier_centers_reduced

  • wannier_centers_reduced (array-like(num_wann,3)) – use the given wannier_centers (reduced) instead of those determined automatically. Incompatible with wannier_centers_cart

  • npar (int) – number of nodes used for parallelization in the __init__ method. Default: multiprocessing.cpu_count()

Notes

  • for tight-binding models it is recommended to use use_wcc_phase = True. In this case the external terms vanish, and

  • one can safely use berry=False, morb=False, and also set ‘external_terms’:False in the parameters of the calculation

set_parameters(**parameters)[source]
set_R_mat(key, value, diag=False, R=None, reset=False, add=False)[source]

Set real-space matrix specified by key. Either diagonal, specific R or full matix. Useful for model calculations

Parameters:
  • key (str) – `SS’, ‘AA’ , etc

  • value (array) –

    • array(num_wann,…) if diag=True . Sets the diagonal part ( if R not set, R=[0,0,0])

    • array(num_wann,num_wann,..) matrix for R (R should be set )

    • array(num_wann,num_wann,nRvec,…)` full spin matrix for all R

    denotes the vector/tensor cartesian dimensions of the matrix element

  • R (list(int)) – list of 3 integer values specifying R. if

  • reset (bool) – allows to reset matrix if it is already set

  • add (bool) – add matrix to the already existing

set_spin(spins, axis=[0, 0, 1], **kwargs)[source]

Set spins along axis in SS(R=0). Useful for model calculations. Note : The spin matrix is purely diagonal, so that <up | sigma_x | down> = 0 For more cversatility use set_R_mat() set_spin_pairs(), set_spin_from_code()

Parameters:
  • spin (one on the following) – 1D array(num_wann) of +1 or -1 spins are along axis

  • axis (array(3)) – spin quantization axis (if spin is a 1D array)

  • **kwargs – optional arguments ‘R’, ‘reset’, ‘add’ see set_R_mat()

set_spin_pairs(pairs)[source]

set SS_R, assuming that each Wannier function is an eigenstate of Sz, :param pairs: list of pair of indices of bands [(up1,down1), (up2,down2), ..] :type pairs: list of tuple :param Notes: :param ——-: :param * For abinitio calculations this is a rough approximation: :param that may be used on own risk.: :param See also set_spin_from_code():

set_spin_from_code(DFT_code='qe')[source]
set SS_R, assuming that each Wannier function is an eigenstate of Sz,

according to the ordering of the ab-initio code

Parameters:
  • DFT_code (str) –

    DFT code used :
    • 'qe' : if bands are grouped by orbital type, in each pair first comes spin-up,then spin-down

    • 'vasp' : if bands are grouped by spin : first come all spin-up, then all spin-down

    1D array(num_wann) of +1 or -1 spins are along axis

  • Notes

  • -------

  • approximation (* This is a rough) –

  • risk (that may be used on own) –

  • Wannier90 (* The pure-spin character may be broken by maximal localization. Recommended to use num_iter=0 in) –

  • name (* if your DFT code has a different) –

  • correspondingly (but uses the same spin ordering as qe or vasp - set DFT_code='qe' or DFT_code='vasp') –

  • ordering (* if your DFT code has a different spin) –

:param use set_spin_pairs():

set_symmetry(symmetry_gen=[])[source]

Set the symmetry group of the System

Parameters:

symmetry_gen (list of symmetry.Symmetry or str) – The generators of the symmetry group.

Notes

  • Only the generators of the symmetry group are essential. However, no problem if more symmetries are provided. The code further evaluates all possible products of symmetry operations, until the full group is restored.

  • Providing Identity is not needed. It is included by default

  • Operations are given as objects of symmetry.Symmetry or by name as str, e.g. 'Inversion' , 'C6z', or products like 'TimeReversal*C2x'.

  • symetyry_gen=[] is equivalent to not calling this function at all

  • Only the point group operations are important. Hence, for non-symmorphic operations, only the rotational part should be given, neglecting the translation.

set_structure(positions, atom_labels, magnetic_moments=None)[source]

Set atomic structure of the system.

Parameters:
  • positions ((num_atom, 3) array_like of float) – Atomic positions in fractional coordinates.

  • atom_labels ((num_atom,) list) – labels (integer, string, etc.) to distinguish species.

  • magnetic_moments ((num_atom, 3) array_like of float (optional)) – Magnetic moment vector of each atom.

set_symmetry_from_structure()[source]

Set the symmetry group of the System. Requires spglib to be installed. System.set_structure() must be called in advance.

For magnetic systems, symmetries involving time reversal are not detected because spglib does not support time reversal symmetry for noncollinear systems.

Symmetrization of the system

System.symmetrize(proj, positions, atom_name, soc=False, magmom=None, DFT_code='qe', method='new')[source]

Symmetrize Wannier matrices in real space: Ham_R, AA_R, BB_R, SS_R,…

Parameters:
  • positions (array) – Positions of each atom.

  • atom_name (list) – Name of each atom.

  • proj (list) –

    Should be the same with projections card in relative Wannier90.win.

    eg: ['Te: s','Te:p']

    If there is hybrid orbital, grouping the other orbitals.

    eg: ['Fe':sp3d2;t2g] Plese don’t use ['Fe':sp3d2;dxz,dyz,dxy]

    ['X':sp;p2] Plese don’t use ['X':sp;pz,py]

  • soc (bool) – Spin orbital coupling.

  • magmom (2D array) – Magnetic momens of each atoms.

  • DFT_code (str) – DFT code used : 'qe' or 'vasp' . This is needed, because vasp and qe have different orbitals arrangement with SOC.(grouped by spin or by orbital type)

  • method (str) – new or old. They give same result but new is faster. old will be eventually removed.

  • Notes – does not update wannier_centers. TODO: make the code update them

From Wannier functions

Wannier90

class wannierberri.system.System_w90(seedname='wannier90', transl_inv=True, guiding_centers=False, fft='fftw', npar=2, kmesh_tol=1e-07, bk_complete_tol=1e-05, **parameters)[source]

Bases: System

System initialized from the Wannier functions generated by Wannier90 code. Reads the .chk, .eig and optionally .mmn, .spn, .uHu, .sIu, and .sHu files

Parameters:
  • seedname (str) – the seedname used in Wannier90

  • transl_inv (bool) – Use Eq.(31) of Marzari&Vanderbilt PRB 56, 12847 (1997) for band-diagonal position matrix elements

  • guiding_centers (bool) – If True, enable overwriting the diagonal elements of the AA_R matrix at R=0 with the Wannier centers calculated from Wannier90.

  • npar (int) – number of processes used in the constructor

  • fft (str) – library used to perform the fast Fourier transform from q to R. fftw or numpy. (practically does not affect performance, anyway mostly time of the constructor is consumed by reading the input files)

  • kmesh_tol (float) – tolerance to consider the b_k vectors (connecting to neighbouring k-points on the grid) belonging to the same shell

  • bk_complete_tol (float) – tolerance to consider the set of b_k shells as complete.

Notes

see also parameters of the System

FPLO

class wannierberri.system.System_fplo(hamdata='+hamdata', **parameters)[source]

Bases: System

System initialized from the +hamdata file written by FPLO code,

Parameters:

hamdata (str) – name (and path) of the “+hamdata” file to be read

Notes

see also parameters of the System

ASE

class wannierberri.system.System_ASE(ase_wannier, ase_R_vectors=False, **parameters)[source]

Bases: System_w90

System initialized from the Wannier functions generated by ASE .

Parameters:

ase_wannier – An object of ASE Wannier .

Notes

see also parameters of the System

From tight-binding models

wannier90_tb.dat file

class wannierberri.system.System_tb(tb_file='wannier90_tb.dat', **parameters)[source]

Bases: System

System initialized from the *_tb.dat file, which can be written either by Wannier90 code, or composed by the user based on some tight-binding model. See Wannier90 code for details of the format.

Parameters:

tb_file (str) – name (and path) of file to be read

Notes

see also parameters of the System

PythTB

class wannierberri.system.System_PythTB(ptb_model, **parameters)[source]

Bases: System_tb_py

This interface is an way to initialize the System class from a tight-binding model created with PythTB. It defines the Hamiltonian matrix Ham_R (from hoppings matrix elements) and the AA_R matrix (from orbital coordinates) used to calculate Berry related quantities.

Parameters:

ptb_model (class) – name of the PythTB tight-binding model class.

Notes

see also parameters of the System_tb_py

TBmodels

class wannierberri.system.System_TBmodels(tbmodel, **parameters)[source]

Bases: System_tb_py

This interface initializes the System class from a tight-binding model created with TBmodels. It defines the Hamiltonian matrix Ham_R (from hoppings matrix elements) and the AA_R matrix (from orbital coordinates) used to calculate Berry related quantities.

Parameters:

tbmodel – name of the TBmodels tight-binding model object.

Notes

see also parameters of the System_tb_py

\(\mathbf{k}\cdot\mathbf{p}\) models

class wannierberri.system.SystemKP(Ham, derHam=None, der2Ham=None, der3Ham=None, kmax=1.0, real_lattice=None, recip_lattice=None, k_vector_cartesian=True, finite_diff_dk=0.0001, maxder=3, **parameters)[source]

Bases: System

A system to describe k.p Hamiltonians. Technically, it is concodered as a periodic system with k-vector limited to the box defined by the reciprocal lattice. a k-vector is always translated to have reduced coordinates between [-1/2,1/2) (In future : translate to the 1BZ for non-simple-cubic lattices)

Parameters:
  • Ham (function) – The Hamiltonian - a funkction of 3D k-vector that returns a (num_waan x num_wann) Hermitean matrix

  • derHam (function) – The cartesian k-derivative of the Hamiltonian - a funkction of 3D k-vector that returns a (num_waan x num_wann x 3) Hermitean (in mn) matrix. If not specified, it will be evaluated numerically from Ham with a finite-difference scheme using the finite_diff_dk parameter.

  • der2Ham (function) – The cartesian second k-derivative of the Hamiltonian - a funkction of 3D k-vector that returns a (num_waan x num_wann x 3 x 3) Hermitean (in mn) matrix If not specified, it will be evaluated numerically from derHam with a finite-difference scheme using the finite_diff_dk parameter.

  • der3Ham (function) – The cartesian second k-derivative of the Hamiltonian - a funkction of 3D k-vector that returns a (num_waan x num_wann x 3 x 3 x 3) Hermitean (in mn) matrix If not specified, it will be evaluated numerically from der2Ham with a finite-difference scheme using the finite_diff_dk parameter.

  • kmax (float) – maximal k-vector (in \(\mathring{\rm A}^{-1}\)) In this case the reciprocal lattice is cubic with size 2kmax

  • real_lattice (array(3,3)) – the lattice vectors of the model (iif kmax is not set)

  • recip_lattice (array(3,3)) – the reciprocal lattice vectors of the model (if ‘kmax’,’real_lattice’ are not set)

  • k_vector_cartesian (bool) – if True, the k-vector in Ham, derHar, der2Ham is given in cartesian coordinates. if False - it is in reduced coordinates (Note : the derivatives are always assumed w.r.t. cartesian coordinates)

  • finite_diff_dk (float) – defines the dk in taking derivatives *in fraction of the reciprocal lattice

Notes

  • if derivatives of hamiltonian are not provided, they are computed with finite diifferences

  • internally, self.Ham and derivatives (self.Ham_cart, self_derHam_cart …) accept k in reduced coordinated.

  • the derivatives are always assumed w.r.t. cartesian coordinates