First - Principles Systems

This specialisation of Z2Pack can handle systems computed with ab initio codes interfacing to Wannier90.

Note

A modified version of Wannier90 is needed for Z2Pack. Please consult the Tutorial for details.

class z2pack.fp.System(*, input_files, kpt_fct, kpt_path, command, executable=None, build_folder='build', file_names=None, mmn_path='wannier90.mmn', num_wcc=None)[source]

Bases: OverlapSystem

System class for systems which are calculated from first principles.

Parameters
  • input_files (list of str) – Paths of the input files.

  • kpt_fct – Function that creates a str specifying the k-points (in the language of the first-principles code used), given a starting_point, last_point, end point and number of k-points N. Can also be a list of functions if k-points need to be written to more than one file.

  • kpt_path (str, or list thereof) – Name of the file where the k-points str belongs. Will append to a file if it matches one of the file_names, and create a separate file else. If kpt_fct is a list, kpt_path should also be a list, specifying the path for each of the functions.

  • command (str) – Command to execute the first principles code.

  • executable (str) – Sets the executable executing the command. If nothing is specified, the subprocess default will be used.

  • build_folder (str) – Folder where the calculation is executed.

  • file_names (list of str) – Names the input files should get in the build_folder. Default behaviour is taking the filenames from the input files.

  • mmn_path (str) – Path to the .mmn output file of Wannier90

  • num_wcc (int) – Number of WCC which should be produced by the system. This parameter can be used to check the consistency of the calculation. By default, no such check is done.

Note

input_files and build_folder can be absolute or relative paths, the rest is relative to build_folder

get_mmn(kpt)[source]

Returns a list of overlap matrices \(M_{m,n}\) corresponding to the given k-points.

Parameters

kpt (list) – The list of k-points for which the overlap matrices are to be computed.

Functions creating k - points input

A collection of functions for creating k-points input for different first-principles codes.

All functions have the same calling structure as prototype().

z2pack.fp.kpoint.abinit(kpt)[source]

Creates a k-point input for ABINIT. It uses kptopt -1 and specifies the k-points string using ndivk and kptbounds.

z2pack.fp.kpoint.prototype(kpt)[source]

Specifies the interface

Parameters

kpt (list of numpy.array) – The list of k-points in the string INCLUDING the final point which should not be in the calculation

z2pack.fp.kpoint.qe(kpt)[source]

Creates a k-point input for Quantum Espresso.

z2pack.fp.kpoint.qe_explicit(kpt)[source]

Creates a k-point input for Quantum Espresso, by explicitly specifying the k-points.

z2pack.fp.kpoint.vasp(kpt)[source]

Creates a k-point input for VASP. It uses the automatic generation scheme with a Gamma centered grid. Note that VASP does not support any kind of k-point line unless they are exactly along one of the reciprocal lattice vectors, and the k-points are evenly spaced.

z2pack.fp.kpoint.wannier90(kpt)[source]

Creates a k-point input for Wannier90. It can be useful when the first-principles code does not generate the k-points in wannier90.win (e.g. with Quantum Espresso).

z2pack.fp.kpoint.wannier90_full(kpt)[source]

Returns both k-point and nearest neighbour input for wannier90.win. This is the recommended function to use for Wannier90 2.1 and higher.

z2pack.fp.kpoint.wannier90_nnkpts(kpt)[source]

Creates the nnkpts input to explicitly specify the nearest neighbours in wannier90.win