QE_methods

libra_py.QE_methods.cryst2cart(a1, a2, a3, r)[source]

Crystal to Cartesian coordinate conversion

An auxiliary function that convert vector <r> in crystal (alat) coordinates with the cell defined by vectors a1, a2, a3, to the Cartesian coordinate xyz

Parameters
  • a1 ([double, double, double]) – one of the unit cell/periodicty vectors

  • a2 ([double, double, double]) – one of the unit cell/periodicty vectors

  • a3 ([double, double, double]) – one of the unit cell/periodicty vectors

  • r ([double, double, double]) – Crystal (fractional) coordinates of the atom

Returns

the Cartesian coordinates of an atom

Return type

[double, double, double]

libra_py.QE_methods.get_QE_normal_modes(filename, verbosity=0)[source]

This function reads the QE phonon calculations output files to get the key information for further normal modes visualization or other types of calculations related to normal modes

Parameters
  • filename (string) – the name of a .dyn file produced by QE code

  • verbosity (int) –

    • 0 - no extra output (default)

    • 1 - print extra stuff

Returns

(Elts, R, U), where:

  • Elts ( list of nat string ): labels of all atoms, nat - is the number of atoms

  • R ( MATRIX(3*nat x 1) ): coordinates of all atoms [Angstrom]

  • U ( MATRIX(ndof x ndof) ): eigenvectors, defining the normal modes

Return type

tuple

Example

>>> get_QE_normal_modes("silicon.dyn1", 1)     # verbose output
>>> get_QE_normal_modes("Cs4SnBr6_T200.dyn1")  # not verbose output
libra_py.QE_methods.out2inp(out_filename, templ_filename, wd, prefix, t0, tmax, dt)[source]

Converts a QE output file with an MD trajectory to a bunch of input files for SCF calculations. These input files all have the same control settings, but differ in atomic coordinates

Parameters
  • out_filename (string) – name of the file which contains the MD trajectory

  • templ_filename (string) – name of the template file for input generation should not contain atomic positions!

  • prefix (string) – the prefix of the files generated as the output

  • wd (string) – working directory where all files will be created/processed - will be created

  • t0 (int) – defines the starting timestep to process (not all the MD timesteps may be precessed)

  • tmax (int) – defines the maximal timestep to process (not all the MD timesteps may be precessed)

  • dt (int) – defines the spacing between frames which are written this is defined as a difference between written configuration indexes. So if dt = 5, the following frames will be written: 0,5,10,15, etc…

Returns

but will create a bunch of new input files in the created working directory

Return type

None

libra_py.QE_methods.out2pdb(out_filename, T, dt, pdb_prefix)[source]

Converts a QE output file with an MD trajectory to a bunch of pdb files

Parameters
  • out_filename (string) – name of the file which contains the MD trajectory. This file is the QE MD trajectory output

  • T (int) – defines the maximal timestep to process (not all the MD timesteps may be precessed)

  • dt (int) – defines the spacing between frames which are written this is defined as a difference between written configuration indexes. So if dt = 5, the following frames will be written: 0,5,10,15, etc…

Returns

but will create a bunch of new pdb files with the trajectory info, including the unit cell params.

Return type

None

Example

>>> QE_methods.out2pdb("x.md.out",250,25,"snaps/snap_")

This will create MD snapshots at times 0 (input configuration), 25 (25-th nuclear configuration), 50, etc. the files will be collcted in the folder /snaps and named snap_0.pdb, snap_25.pdb, snap_50.pdb, etc.

libra_py.QE_methods.out2xyz(out_filename, T, dt, xyz_filename)[source]

This function converts the QE output file into a .xyz trajectory file. No more than T steps from the out_filename file will be used

Parameters
  • out_filename (string) – name of the file which contains the MD trajectory. This file is the QE MD trajectory output

  • T (int) – defines the maximal timestep to process (not all the MD timesteps may be precessed)

  • dt (int) – defines the spacing between frames which are written this is defined as a difference between written configuration indexes. So if dt = 5, the following frames will be written: 0,5,10,15, etc…

  • xyz_filename (string) – the name of the file to which the resulting .xyz trajectory will be written

Returns

but will create a .xyz file with the trajectory.

Return type

None

Example

>>> QE_methods.out2xyz("x.md.out",250,25,"snaps/traj.xyz")

This will create the MD trajectory file in .xyz format with the snapshots takes at times 0 (input configuration), 25 (25-th nuclear configuration), 50, etc. the snapshots will written in the file trahj.xyz in the folder /snaps

libra_py.QE_methods.read_all(params)[source]

This function reads index, wfc and grid files from a given directory The number of wfc and grid files may be larger than 1 - this is the case of spin-polarized or multiple k-points calculations

Args:

params ( dictionary ): Parameters controlling the simulation parameters

  • params[“wd”] ( string ): the name of a “working directory (can be removed once the calculatons

    are done)” that will be created during this function execution - this is where the temporary files are written to [default: “wd”]

  • params[“prefix”] ( string ): the location of the folder containing index.xml, wfc.*, and grid.* files [default: “x0.export” ]

  • params[“read_wfc”] ( 0 or 1 ): whether or not to read the wfc coefficients. [ default: 1 ]

  • params[“read_grid”] ( 0 or 1 ): whether or not to read the grid informations. [ default: 1 ]

  • params[“verb0”] ( 0 or 1 ): turn off/on the extra printout while reading index.xml. [ default: 0 ]

  • params[“verb1”] ( 0 or 1 ): turn off/on the extra printout while reading wfc.*. [ default: 0 ]

  • params[“verb2”] ( 0 or 1 ): turn off/on the extra printout while reading grid.*. [ default: 0 ]

  • params[“nac_method”] ( 0, 1, 2 ): the expectations about what format to read:

    • 0 - non-SOC, non-polarized

    • 1 - non-SOC, spin-polarized

    • 2 - SOC, non-collinear

  • params[“minband”] ( int ): index of the lowest energy orbital to include

    in the active space, counting starts from 1 [ default: 1]

  • params[“maxband”] ( int ): index of the highest energy orbital to include

    in the active space, counting starts from 1 [ defaults: 2]

Returns

( info, e, coeff, grid ), where

  • info ( dictionary ): general descritor info ..seealso::`QE_methods.read_qe_index`

  • e ( list of CMATRIX(norbs, norbs) ): band energies for each k-pints ..seealso::`QE_methods.read_qe_index`

  • coeff ( list of CMATRIX(npw, len(act_space)) objects ): such the

    coeff[k] are the MOs in the plane wave basis for the k-point k

  • grid ( list of VECTOR objects ): the grid point vectors [ units: tpiba ]

The number of elements in each list is determined by the number of k points Note that, for spin-polarized calculations, the number of k-points is always twice that of the non-spin-polarized or non-collinear k-points

Return type

tuple

libra_py.QE_methods.read_info(params)[source]

This fucntions reads the output from QE calculations, and stores the output information in dictionaries

Parameters

params (dictionary) –

Calculation control parameters

  • params[“wd”] ( string ): the name of a “working directory (can be removed once the calculatons

    are done)” that will be created during this function execution - this is where the temporary files are written to [default: “wd”]

  • params[“nac_method”] ( int ): selects the type of output to analyze:

    • 0 : non-spin-polarized calculations

    • 1 : spin-polarized calculations

    • 2 : non-collinear calculation (SOC) only

    • 3 : spin-polarized and non-collinear calculation (SOC)

Returns

( info0, all_e_dum0, info1, all_e_dum1 ):

info0 ( dictionary ): QE calculations info for the spin-diabatic calculations

all_e_dum0 ( list of CMATRIX(norb, norb) objects ): (eigen)energies for all the k-points for

the spin-diabatic calculations

info1 ( dictionary ): QE calculations info for the non-collinear (spin-adiabatic) calculations

all_e_dum1 ( list of CMATRIX(norb, norb) objects ): (eigen)energies for all the k-points for

the non-collinear (spin-adiabatic) calculations

..seealso:: `QE_methods.read_qe_index`

Return type

tuple

libra_py.QE_methods.read_md_data(filename)[source]

Read in the QE MD data stored in an XML file

Parameters

filename (string) – the name of the xml file that contains an MD data this function is specifically tailored for the QE output format

Returns

(R, V, A, M, E), where:

  • R ( MATRIX(ndof x nsteps-1) ): coordinates of all DOFs for all mid-timesteps [Bohr]

  • V ( MATRIX(ndof x nsteps-1) ): velocities of all DOFs for all mid-timesteps [a.u. of velocity]

  • A ( MATRIX(ndof x nsteps-1) ): accelerations of all DOFs for all mid-timesteps [a.u. of acceleration]

  • M ( MATRIX(ndof x 1) ): masses of all DOFs [a.u. of mass]

  • E (list of ndof/3): atom names (elements) of all atoms

Return type

tuple

libra_py.QE_methods.read_md_data_cell(filename)[source]

Read in the QE MD unit cell vectors stored in an XML file

Parameters

filename (string) – the name of the xml file that contains an MD data this function is specifically tailored for the QE output format

Returns

cell coordinates for all timesteps, in the format [Bohr]:

C(0,0) = a1.x C(0,1) = a1.y C(0,2) = a1.z C(1,0) = a2.x C(1,1) = a2.y C(1,2) = a2.z C(2,0) = a3.x C(2,1) = a3.y C(2,2) = a3.z

Return type

MATRIX(9 x nsteps)

libra_py.QE_methods.read_md_data_xyz(filename, PT, dt)[source]

Read in the MD trajectory stored in an XYZ format

Parameters
  • filename (string) – the name of the xyz file that contains an MD data

  • ( dict{ string (PT) – float } ): definition of the masses of different atomic species [masses in amu, where 1 is the mass of H]

  • dt (double) – MD nuclear integration time step [a.u. of time]

Returns

(R, V, M, E), where:

  • R ( MATRIX(ndof x nsteps-1) ): coordinates of all DOFs for all mid-timesteps [Bohr]

  • V ( MATRIX(ndof x nsteps-1) ): velocities of all DOFs for all mid-timesteps [a.u. of velocity]

  • M ( MATRIX(ndof x 1) ): masses of all DOFs [a.u. of mass]

  • E (list of ndof/3): atom names (elements) of all atoms

Return type

tuple

libra_py.QE_methods.read_md_data_xyz2(filename, PT)[source]

Read in the MD trajectory stored in an XYZ format This version does not compute velocities - only gets coordinates

Parameters
  • filename (string) – the name of the xyz file that contains an MD data

  • ( dict{ string (PT) – float } ): definition of the masses of different atomic species [masses in amu, where 1 is the mass of H]

Returns

(R, E), where:

  • R ( MATRIX(ndof x nsteps-1) ): coordinates of all DOFs for all mid-timesteps [Bohr]

  • E (list of ndof/3): atom names (elements) of all atoms

Return type

tuple

libra_py.QE_methods.read_qe_index(filename, orb_list, verbose=0)[source]

This functions reads an ASCII/XML format file containing wavefunction and returns the coefficients of the plane waves that constitute the wavefunction

Parameters
  • filename (string) – This is the name of the file we will be reading to construct a wavefunction

  • orb_list (list of ints) – The indices of the orbitals which we want to consider. Orbitals indexing at 1, not 0

  • verbose (int) –

    The flag controlling the amout of extra output:

    • 0 - no extra output (default)

    • 1 - print extra stuff

Returns

( info, all_e ), where:

  • info ( dictionary ): contains all the descritive information about the QE calculations and the system

    It contains the following info:

    • info[“nspin”] ( int ): the type of calculations: 1 - non-polarized, 2 - polarized, 4 - non-collinear

    • info[“nk”] ( int ): the number of k-points

    • info[“nbnd”] ( int ): the number of orbitals in each k-point

    • info[“efermi”] ( double ): the Fermi energy [units: a.u.]

    • info[“alat”] ( double ): lattice constant, a [units: Bohr]

    • info[“omega”] ( double ): unit cell volume [units: Bohr^3]

    • info[“tpiba”] ( double ): reciprocal cell size, 2*pi/a [units: Bohr^-1]

    • info[“tpiba2”] ( double ): squared reciprocal cell size, (2*pi/a)^2 [units: Bohr^-2]

    • info[“a1”] ( VECTOR ): direct lattice vector 1 [units: Bohr]

    • info[“a2”] ( VECTOR ): direct lattice vector 2 [units: Bohr]

    • info[“a3”] ( VECTOR ): direct lattice vector 3 [units: Bohr]

    • info[“b1”] ( VECTOR ): reciprocal lattice vector 1 [units: Bohr^-1]

    • info[“b2”] ( VECTOR ): reciprocal lattice vector 2 [units: Bohr^-1]

    • info[“b3”] ( VECTOR ): reciprocal lattice vector 3 [units: Bohr^-1]

    • info[“weights”] ( list ): weights of all the k-points in evaluating the integrals

    • info[“k”] ( list of VECTOR objects ): coordinates of the k-points [units: tpiba]

all_e ( list of CMATRIX(norb, norb) objects ): orbital energies for all the k-points, such that

all_e[k] is a diagonal matrix (complex-valued) of orbital energies - only for those that are of interest Here, norb = len(orb_list) and the matrix will only contain numbers that correspond to orbitals defined in the orb_list argument [units: a.u.]

Return type

tuple

libra_py.QE_methods.read_qe_schema(filename, verbose=0)[source]

This functions reads an ASCII/XML format file containing basic info about the QE run

Parameters
  • filename (string) – This is the name of the file we will be read [ usually: “x0.save/data-file-schema.xml”]

  • verbose (int) –

    The flag controlling the amout of extra output:

    • 0 - no extra output (default)

    • 1 - print extra stuff

Returns

( info ): the dictionary containaining the descritive information about the QE calculations and the system

It contains the following info:

  • info[“conv”] ( int ): the number of steps to converge SCF, 0 - means no convergence

  • info[“etot”] ( double ): the total energy (electronic + nuclear) [ units: Ha ]

  • info[“nbnd”] ( int ): the number of bands

  • info[“nelec”] ( int ): the number of electrons

  • info[“efermi”] ( double ): the Fermi energy [ units: Ha ]

  • info[“alat”] ( double ): lattice constant [ units: Bohr ]

  • info[“nat”] ( int ): the number of atoms

  • info[“coords”] ( MATRIX(3*nat, 1) ): the atomic coordinates [ units: Bohr ]

  • info[“atom_labels”] ( list of nat strings ): the atomic names

  • info[“forces”] ( MATRIX(3*nat, 1) ): the atomic forces [ units: Ha/Bohr ]

Return type

dictionary

libra_py.QE_methods.read_qe_wfc(filename, orb_list, verbose=0)[source]

This functions reads an ASCII/XML format file containing wavefunction and returns the coefficients of the plane waves that constitute the wavefunction

Parameters
  • filename (string) – This is the name of the file we will be reading to construct a wavefunction

  • orb_list (list of ints) – The indices of the orbitals which we want to consider. Orbitals indexing at 1, not 0

  • verbose (int) –

    The flag controlling the amout of extra output:

    • 0 - no extra output (default)

    • 1 - print extra stuff

Returns

The plane-wave expansion coefficients for all orbitals

Return type

CMATRIX(ngw,norbs)

libra_py.QE_methods.read_qe_wfc_grid(filename, verbose=0)[source]

This functions reads an ASCII/XML format file containing grid of G-points for given k-point and returns a list of VECTOR objects

Parameters
  • filename (string) – This is the name of the file we will be reading to construct a wavefunction

  • verbose (int) –

    The flag controlling the amout of extra output:

    • 0 - no extra output (default)

    • 1 - print extra stuff

Returns

the definitions of all G points in the present calculation

Return type

VectorList = list of VECTOR objects

libra_py.QE_methods.read_qe_wfc_info(filename, verbose=0)[source]

This functions reads an ASCII/XML format file containing wavefunction and returns some descriptors

Parameters
  • filename (string) – This is the name of the file we will be reading to construct a wavefunction

  • verbose (int) –

    The flag controlling the amout of extra output:

    • 0 - no extra output (default)

    • 1 - print extra stuff

Returns

a dictionary object that will contain key-value pairs with the descriptors

  • res[“ngw”] ( int ): the number of plane waves needed to represent the orbital

  • res[“igwx”] ( int ): the number of the G points = plane waves needed to

    represent the orbital for given k-point. Use this number when working with multiple k-points

  • res[“nbnd”] ( int ): the number of bands (orbitals)

  • res[“nspin”] ( int ): the desriptor of the spin-polarisation in the calculation

    • 1: unpolarized

    • 2: polarized

    • 4: non-collinear

  • res[“gamma_only”]: ( string = “T” or “F” ): T - use the Gamma-point storae trick, F otherwise

  • res[“ik”] ( int ): index of the k point wfc

  • res[“nk”] ( int ): the number of k-points in the wfc

Return type

dictionary

libra_py.QE_methods.read_wfc_grid(params)[source]

Read the coefficients and energies for the multi k-points cases, even if some cases require gamma only

Parameters

params (dictionary) –

The control parameters, which may contain:

  • param[“nac_method”] ( 0, 1, 2, 3 ): the expectations about what format to read:

    • 0 : non-spin-polarized calculations [ default ]

    • 1 : spin-polarized calculations

    • 2 : non-collinear calculation (SOC) only

    • 3 : spin-polarized and non-collinear calculation (SOC)

  • params[“minband”] ( int ): index of the lowest energy orbital to include

    in the active space in the non-SOC (spin-diabatic) calculations, counting starts from 1. Used for nac_method == 0, 1, 3. [ default: 1]

  • params[“maxband”] ( int ): index of the highest energy orbital to include

    in the active space in the non-SOC (spin-diabatic) calculations, counting starts from 1. Used for nac_method == 0, 1, 3. [ defaults: 2]

  • params[“minband_soc”] ( int ): index of the lowest energy orbital pair to include

    in the active space in the SOC (spin-adiabatic) calculations, counting starts from 1. Used for nac_method == 2 and 3. [ default: 1]

  • params[“maxband_soc”] ( int ): index of the highest energy orbital pair to include

    in the active space in the SOC (spin-adiabatic) calculations, counting starts from 1. Used for nac_method == 2 and 3. [ defaults: 2]

Returns

( res_curr, res_next ), Here _curr, refers to the current timestep properties,

and the _next refers to the consecutive timestep properties. Each element of the output is a dictionary with the following elements:

  • res_curr[“Coeff_dia”] ( list of CMATRIX(npw, len(act_space)) objects )the

    wavefunction coefficients in the planewave basis for the spin-diabatic wavefunctions, such that res_curr[“Coeff_dia”][k] is a matrix for the k-point with index k. Only for nac_method == 0, 1, and 3.

  • res_curr[“E_dia”] ( list of MATRIX(len(act_space), len(act_space)) objects )the MO

    energies for the spin-diabatic wavefunctions, such that res_curr[“E_dia”][k] is a matrix for the k-point with index k. Only for nac_method == 0, 1, and 3.

  • res_curr[“Coeff_adi”] ( list of CMATRIX(npw, len(act_space)) objects )the

    wavefunction coefficients in the planewave basis for the spin-adiabatic wavefunctions, such that res_curr[“Coeff_adi”][k] is a matrix for the k-point with index k. Only for nac_method == 2 and 3.

  • res_curr[“E_adi”] ( list of MATRIX(len(act_space), len(act_space)) objects )the MO

    energies for the spin-adiabatic wavefunctions, such that res_curr[“E_adi”][k] is a matrix for the k-point with index k. Only for nac_method == 2 and 3.

  • res_curr[“grid”] ( list of VECTOR objects ): the grid point vectors [ units: tpiba ]

Return type

tuple

libra_py.QE_methods.run_qe(params, t, dirname0, dirname1)[source]

This function runs necessary QE calculations as defined by the “params” dictionary

Parameters
  • params (dictionary) –

    A dictionary containing important simulation parameters

    • params[“BATCH_SYSTEM”] ( string ): the name of the job submission command

      use “srun” if run calculations on SLURM system or “mpirun” if run on PBS system [default: “srun”]

    • params[“NP”] ( int ): the number of nodes on which execute calculations

      [default: 1]

    • params[“EXE”] ( string ): the name of the program to be executed. This may be

      the absolute path to the QE (pw.x) binary

    • params[“EXE_EXPORT”] ( string ): the name of the program that converts the binary files

      with the QE wavefunctions to the text format (pw_export.x). The name includes the absolute path to the binary

    • params[“prefix0”] ( string ): the name of scf template input file - it should

      contain all the parameters controlling the computational methodology for QE. If the file is called “x0.scf.in”, use “x0.scf” as the value of the “prefix0” [default: “x0.scf”]

    • params[“prefix1”] ( string ): the name of scf template input file - it should

      contain all the parameters controlling the computational methodology for QE. Presently is used for SOC-enabled calculations, whereas the “prefix0” defines the no-SOC calculations. If the file is called “x1.scf.in”, use “x1.scf” as the value of the “prefix1” [default: “x1.scf”]

    • params[“nac_method”] ( int ): selects the type of calculations to perform:

      • 0 : non-spin-polarized calculations (needs only “prefix0”)

      • 1 : spin-polarized calculations (needs only “prefix0”)

      • 2 : non-collinear calculation (SOC) only (needs only “prefix1”)

      • 3 : spin-polarized and non-collinear calculation (SOC) (needs both “prefix0” and “prefix1”)

      [default: 0]

    • params[“wd”] ( string ): the name of a “working directory (can be removed once the calculatons

      are done)” that will be created during this function execution.

  • t (int) – the current time step

  • dirname0 (string) – Name of the temporary directory where data will be stored for the case without SOC

  • dirname1 (string) – Name of the temporary directory where data will be stored for the case with SOC

libra_py.QE_methods.xyz2inp(out_filename, templ_filename, wd, prefix, t0, tmax, dt)[source]

Converts a xyz trajectory file with an MD trajectory to a bunch of input files for SCF calculations. These input files all have the same control settings, but differ in atomic coordinates

Parameters
  • out_filename (string) – name of the file which contains the MD trajectory (in xyz format)

  • templ_filename (string) – name of the template file for input generation should not contain atomic positions!

  • wd (string) – working directory where all files will be created/processed - will be created

  • prefix (string) – the prefix of the files generated as the output

  • t0 (int) – defines the starting timestep to process (not all the MD timesteps may be precessed)

  • tmax (int) – defines the maximal timestep to process (not all the MD timesteps may be precessed)

  • dt (int) – defines the spacing between frames which are written this is defined as a difference between written configuration indexes. So if dt = 5, the following frames will be written: 0,5,10,15, etc…

Returns

but will create a bunch of new input files in the created working directory

Return type

None