pymatgen.io.lobster.outputs module
Module for reading Lobster output files. For more information on LOBSTER see www.cohp.de.
- class Bandoverlaps(filename: str = 'bandOverlaps.lobster')[source]
Bases:
object
Class to read in bandOverlaps.lobster files. These files are not created during every Lobster run. .. attribute: bandoverlapsdict is a dict of the following form:
{spin:{“kpoint as string”: {“maxDeviation”: float that describes the max deviation, “matrix”: 2D array of the size number of bands times number of bands including the overlap matrices with } }}
- Parameters:
filename – filename of the “bandOverlaps.lobster” file
- has_good_quality_check_occupied_bands(number_occ_bands_spin_up: int, number_occ_bands_spin_down: int | None = None, spin_polarized: bool = False, limit_deviation: float = 0.1) bool [source]
Will check if the deviation from the ideal bandoverlap of all occupied bands is smaller or equal to limit_deviation
Args: number_occ_bands_spin_up (int): number of occupied bands of spin up number_occ_bands_spin_down (int): number of occupied bands of spin down spin_polarized (bool): If True, then it was a spin polarized calculation limit_deviation (float): limit of the maxDeviation :returns: Boolean that will give you information about the quality of the projection
- has_good_quality_maxDeviation(limit_maxDeviation: float = 0.1) bool [source]
Will check if the maxDeviation from the ideal bandoverlap is smaller or equal to limit_maxDeviation :param limit_maxDeviation: limit of the maxDeviation
- Returns:
Boolean that will give you information about the quality of the projection
- class Charge(filename: str = 'CHARGE.lobster')[source]
Bases:
object
Class to read CHARGE files generated by LOBSTER
- Parameters:
filename – filename for the CHARGE file, typically “CHARGE.lobster”
- class Cohpcar(are_coops: bool = False, are_cobis: bool = False, filename: str | None = None)[source]
Bases:
object
Class to read COHPCAR/COOPCAR files generated by LOBSTER.
- Parameters:
are_coops – Determines if the file is a list of COHPs or COOPs. Default is False for COHPs.
are_cobis – Determines if the file is a list of COHPs or COOPs. Default is False for COHPs.
filename – Name of the COHPCAR file. If it is None, the default file name will be chosen, depending on the value of are_coops.
- class Doscar(doscar: str = 'DOSCAR.lobster', structure_file: str = 'POSCAR', dftprogram: str = 'Vasp')[source]
Bases:
object
Class to deal with Lobster’s projected DOS and local projected DOS. The beforehand quantum-chemical calculation was performed with VASP
- pdos[source]
- List of Dict including numpy arrays with pdos. Access as pdos[atomindex]['orbitalstring']['Spin.up/Spin.down']
- energies[source]
- numpy array of the energies at which the DOS was calculated (in eV, relative to Efermi)
- tdensities[source]
- tdensities[Spin.up]: numpy array of the total density of states for the Spin.up contribution at each of the
- energies
- tdensities[Spin.down]: numpy array of the total density of states for the Spin.down contribution at each of the
- energies
if is_spin_polarized=False: tdensities[Spin.up]: numpy array of the total density of states
- itdensities:
- itdensities[Spin.up]: numpy array of the total density of states for the Spin.up contribution at each of the
- energies[source]
- itdensities[Spin.down]: numpy array of the total density of states for the Spin.down contribution at each of the
- energies
if is_spin_polarized=False: itdensities[Spin.up]: numpy array of the total density of states
- Parameters:
doscar – DOSCAR filename, typically “DOSCAR.lobster”
structure_file – for vasp, this is typically “POSCAR”
dftprogram – so far only “vasp” is implemented
- property completedos: LobsterCompleteDos[source]
CompleteDos
- Type:
return
- class Fatband(filenames='.', vasprun='vasprun.xml', Kpointsfile='KPOINTS')[source]
Bases:
object
Reads in FATBAND_x_y.lobster files
- Parameters:
filenames (list or string) – can be a list of file names or a path to a folder folder from which all “FATBAND_*” files will be read
vasprun – corresponding vasprun file
Kpointsfile – KPOINTS file for bandstructure calculation, typically “KPOINTS”
- class Grosspop(filename: str = 'GROSSPOP.lobster')[source]
Bases:
object
Class to read in GROSSPOP.lobster files.
- Parameters:
filename – filename of the “GROSSPOP.lobster” file
- get_structure_with_total_grosspop(structure_filename: str) Structure [source]
Get a Structure with Mulliken and Loewdin total grosspopulations as site properties :param structure_filename: filename of POSCAR :type structure_filename: str
- Returns:
Structure Object with Mulliken and Loewdin total grosspopulations as site properties
- class Icohplist(are_coops: bool = False, are_cobis: bool = False, filename: str | None = None)[source]
Bases:
object
Class to read ICOHPLIST/ICOOPLIST files generated by LOBSTER.
- Parameters:
are_coops – Determines if the file is a list of ICOOPs. Defaults to False for ICOHPs.
are_cobis – Determines if the file is a list of ICOBIs. Defaults to False for ICOHPs.
filename – Name of the ICOHPLIST file. If it is None, the default file name will be chosen, depending on the value of are_coops.
- class Lobsterout(filename='lobsterout')[source]
Bases:
object
Class to read in the lobsterout and evaluate the spilling, save the basis, save warnings, save infos
- Parameters:
filename – filename of lobsterout
- class MadelungEnergies(filename: str = 'MadelungEnergies.lobster')[source]
Bases:
object
Class to read MadelungEnergies.lobster files generated by LOBSTER
- Parameters:
filename – filename of the “MadelungEnergies.lobster” file
- class SitePotential(filename: str = 'SitePotentials.lobster')[source]
Bases:
object
Class to read SitePotentials.lobster files generated by LOBSTER
- Parameters:
filename – filename for the SitePotentials file, typically “SitePotentials.lobster”
- class Wavefunction(filename, structure)[source]
Bases:
object
Class to read in wave function files from Lobster and transfer them into an object of the type VolumetricData
- Parameters:
filename – filename of wavecar file from Lobster
structure – Structure object (e.g., created by Structure.from_file(“”))
- get_volumetricdata_density()[source]
Will return a VolumetricData object including the imaginary part of the wave function
Returns: VolumetricData object
- get_volumetricdata_imaginary()[source]
Will return a VolumetricData object including the imaginary part of the wave function
Returns: VolumetricData object
- get_volumetricdata_real()[source]
Will return a VolumetricData object including the real part of the wave function
Returns: VolumetricData object
- set_volumetric_data(grid, structure)[source]
Will create the VolumetricData Objects
- Parameters:
grid – grid on which wavefunction was calculated, e.g. [1,2,2]
structure – Structure object
- write_file(filename='WAVECAR.vasp', part='real')[source]
Will save the wavefunction in a file format that can be read by VESTA This will only work if the wavefunction from lobster was constructed with: “printLCAORealSpaceWavefunction kpoint 1 coordinates 0.0 0.0 0.0 coordinates 1.0 1.0 1.0 box bandlist 1 2 3 4 5 6 ” or similar (the whole unit cell has to be covered!)
- Parameters:
filename – Filename for the output, e.g., WAVECAR.vasp
part – which part of the wavefunction will be saved (“real” or “imaginary”)