pymatgen.io.lobster.inputs module¶
Module for reading Lobster output files. For more information on LOBSTER see www.cohp.de.
-
class
Lobsterin
(settingsdict: dict)[source]¶ Bases:
dict
,monty.json.MSONable
This class can handle and generate lobsterin files Furthermore, it can also modify INCAR files for lobster, generate KPOINT files for fatband calculations in Lobster, and generate the standard primitive cells in a POSCAR file that are needed for the fatband calculations. There are also several standard lobsterin files that can be easily generated.
- Parameters
settingsdict – dict to initialize Lobsterin
-
AVAILABLEKEYWORDS
= ['COHPstartEnergy', 'COHPendEnergy', 'basisSet', 'cohpGenerator', 'gaussianSmearingWidth', 'saveProjectionToFile', 'basisfunctions', 'skipdos', 'skipcohp', 'skipcoop', 'skipPopulationAnalysis', 'skipGrossPopulation', 'userecommendedbasisfunctions', 'loadProjectionFromFile', 'forceEnergyRange', 'DensityOfEnergy', 'BWDF', 'BWDFCOHP', 'skipProjection', 'createFatband', 'writeBasisFunctions', 'writeMatricesToFile', 'realspaceHamiltonian', 'realspaceOverlap', 'printPAWRealSpaceWavefunction', 'printLCAORealSpaceWavefunction', 'noFFTforVisualization', 'RMSp', 'onlyReadVasprun.xml', 'noMemoryMappedFiles', 'skipPAWOrthonormalityTest', 'doNotIgnoreExcessiveBands', 'doNotUseAbsoluteSpilling', 'skipReOrthonormalization', 'forceV1HMatrix', 'useOriginalTetrahedronMethod', 'useDecimalPlaces', 'kSpaceCOHP'][source]¶
-
BOOLEANKEYWORDS
= ['saveProjectionToFile', 'skipdos', 'skipcohp', 'skipcoop', 'loadProjectionFromFile', 'forceEnergyRange', 'DensityOfEnergy', 'BWDF', 'BWDFCOHP', 'skipPopulationAnalysis', 'skipGrossPopulation', 'userecommendedbasisfunctions', 'skipProjection', 'writeBasisFunctions', 'writeMatricesToFile', 'noFFTforVisualization', 'RMSp', 'onlyReadVasprun.xml', 'noMemoryMappedFiles', 'skipPAWOrthonormalityTest', 'doNotIgnoreExcessiveBands', 'doNotUseAbsoluteSpilling', 'skipReOrthonormalization', 'forceV1HMatrix', 'useOriginalTetrahedronMethod', 'forceEnergyRange', 'bandwiseSpilling', 'kpointwiseSpilling'][source]¶
-
FLOATKEYWORDS
= ['COHPstartEnergy', 'COHPendEnergy', 'gaussianSmearingWidth', 'useDecimalPlaces', 'COHPSteps'][source]¶
-
STRINGKEYWORDS
= ['basisSet', 'cohpGenerator', 'realspaceHamiltonian', 'realspaceOverlap', 'printPAWRealSpaceWavefunction', 'printLCAORealSpaceWavefunction', 'kSpaceCOHP'][source]¶
-
diff
(other)[source]¶ Diff function for lobsterin. Compares two lobsterin and indicates which parameters are the same. Similar to the diff in INCAR. :param other: Lobsterin object to compare to :type other: Lobsterin
- Returns
dict with differences and similarities
-
classmethod
from_file
(lobsterin: str)[source]¶ - Parameters
lobsterin (str) – path to lobsterin
- Returns
Lobsterin object
-
static
get_all_possible_basis_functions
(structure: pymatgen.core.structure.Structure, potcar_symbols: list, address_basis_file_min: str = '/Users/shyue/repos/pymatgen/pymatgen/io/lobster/lobster_basis/BASIS_PBE_54_min.yaml', address_basis_file_max: str = '/Users/shyue/repos/pymatgen/pymatgen/io/lobster/lobster_basis/BASIS_PBE_54_max.yaml')[source]¶ - Parameters
structure – Structure object
potcar_symbols – list of the potcar symbols
address_basis_file_min – path to file with the minium required basis by the POTCAR
address_basis_file_max – path to file with the largest possible basis of the POTCAR
Returns: List of dictionaries that can be used to create new Lobsterin objects in standard_calculations_from_vasp_files as dict_for_basis
-
static
get_basis
(structure: pymatgen.core.structure.Structure, potcar_symbols: list, address_basis_file: str = '/Users/shyue/repos/pymatgen/pymatgen/io/lobster/lobster_basis/BASIS_PBE_54_standard.yaml')[source]¶ will get the basis from given potcar_symbols (e.g., [“Fe_pv”,”Si”] #include this in lobsterin class :param structure: Structure object :type structure: Structure :param potcar_symbols: list of potcar symbols
- Returns
returns basis
-
classmethod
standard_calculations_from_vasp_files
(POSCAR_input: str = 'POSCAR', INCAR_input: str = 'INCAR', POTCAR_input: Optional[str] = None, dict_for_basis: Optional[dict] = None, option: str = 'standard')[source]¶ will generate Lobsterin with standard settings
- Parameters
POSCAR_input (str) – path to POSCAR
INCAR_input (str) – path to INCAR
POTCAR_input (str) – path to POTCAR
dict_for_basis (dict) – can be provided: it should look the following: dict_for_basis={“Fe”:’3p 3d 4s 4f’, “C”: ‘2s 2p’} and will overwrite all settings from POTCAR_input
option (str) – ‘standard’ will start a normal lobster run where COHPs, COOPs, DOS, CHARGE etc. will be calculated ‘standard_from_projection’ will start a normal lobster run from a projection ‘standard_with_fatband’ will do a fatband calculation, run over all orbitals ‘onlyprojection’ will only do a projection ‘onlydos’ will only calculate a projected dos ‘onlycohp’ will only calculate cohp ‘onlycoop’ will only calculate coop ‘onlycohpcoop’ will only calculate cohp and coop
- Returns
Lobsterin Object with standard settings
-
write_INCAR
(incar_input: str = 'INCAR', incar_output: str = 'INCAR.lobster', poscar_input: str = 'POSCAR', isym: int = - 1, further_settings: Optional[dict] = None)[source]¶ Will only make the run static, insert nbands, make ISYM=-1, set LWAVE=True and write a new INCAR. You have to check for the rest. :param incar_input: path to input INCAR :type incar_input: str :param incar_output: path to output INCAR :type incar_output: str :param poscar_input: path to input POSCAR :type poscar_input: str :param isym: isym equal to -1 or 0 are possible. Current Lobster version only allow -1. :type isym: int :param further_settings: A dict can be used to include further settings, e.g. {“ISMEAR”:-5} :type further_settings: dict
-
static
write_KPOINTS
(POSCAR_input: str = 'POSCAR', KPOINTS_output='KPOINTS.lobster', reciprocal_density: int = 100, isym: int = - 1, from_grid: bool = False, input_grid: list = [5, 5, 5], line_mode: bool = True, kpoints_line_density: int = 20, symprec: float = 0.01)[source]¶ writes a KPOINT file for lobster (only ISYM=-1 and ISYM=0 are possible), grids are gamma centered :param POSCAR_input: path to POSCAR :type POSCAR_input: str :param KPOINTS_output: path to output KPOINTS :type KPOINTS_output: str :param reciprocal_density: Grid density :type reciprocal_density: int :param isym: either -1 or 0. Current Lobster versions only allow -1. :type isym: int :param from_grid: If True KPOINTS will be generated with the help of a grid given in input_grid. Otherwise,
they will be generated from the reciprocal_density
- Parameters
input_grid (list) – grid to generate the KPOINTS file
line_mode (bool) – If True, band structure will be generated
kpoints_line_density (int) – density of the lines in the band structure
symprec (float) – precision to determine symmetry
-
static
write_POSCAR_with_standard_primitive
(POSCAR_input='POSCAR', POSCAR_output='POSCAR.lobster', symprec=0.01)[source]¶ writes a POSCAR with the standard primitive cell. This is needed to arrive at the correct kpath :param POSCAR_input: filename of input POSCAR :type POSCAR_input: str :param POSCAR_output: filename of output POSCAR :type POSCAR_output: str :param symprec: precision to find symmetry :type symprec: float