overlap_population

The write_lobsterin function

write_lobsterin(directory='.', adsorbate_atoms=['C', 'H', 'O', 'N'], basisSet='pbeVaspFit2015', basisfunctions={'C': '2s 2p', 'H': '1s', 'N': '2s 2p', 'O': '2s 2p', 'Pt': '5d 6s'}, gaussianSmearingWidth=0.003)[source]

Write a lobster input file

Parameters
  • directory (str) – directory to write the lobsterin file. Must also contain a CONTCAR and DOSCAR file

  • adsorbate_atoms (list[str or int]) – adsorbate atom symbols or indices

  • basisSet (str) – basis used to projected density onto orbitals

  • basisfunctions (dict) – dictionary of atomic symbols and corresponding basis functions

  • gaussianSmearingWidth (float) – float if Gaussian smearing is used. If tetrahedron method is used then set to None

The get_example_data function

get_example_data()[source]

Get default path to experimental crystal ovelap populaiton data

Returns

data_path – path to example lobster data

Return type

str

The get_all_lobster_files function

get_all_lobster_files(directory, file_type='COOPCAR.lobster')[source]

Get all DOSCAR and CONTCAR file paths

Parameters
  • directory (str) – root directory to look for DOSCAR files

  • file_type (string) – file type for which to search and return file path

Returns

lobster_files – list of paths to lobster files of type file_type

Return type

list[str]

The get_bonding_states function

get_bonding_states(orbital_indices, dos_energies, pcoop, pcoop_energies, set_antibonding_zero=False, emax=inf)[source]

method for obtaining bonding fraction of dos or dos-like array

Parameters
  • orbital_indices (list[list]) – list of orbital indices for each molecular orbital

  • dos_energies (numpy.ndarray) – energies at which dos is calculated

  • pcoop (numpy.ndarray) – projected orbital overlap populations

  • pcoop_energies (numpy.ndarray) – energies at which the pcoop is calculated

  • set_antibonding_zero (bool) – if true, set antibonding populations to zero to look at total instead of net bonding characteristics

Returns

dos_bonding – 1-D array of the relative bonding for either dos-like array

Return type

numpy.ndarray

The OVERLAP_POPULATION class

class OVERLAP_POPULATION(file_name='COOPCAR.lobster')[source]

Class for analyzing overlap population bonding data

Parameters

file_name (str) – full lobster file location

Variables
  • emax (float) – maximum energy level

  • emin (float) – minimum energy level

  • ndos (int) – number of descritized energy levels

  • e_fermi (float) – highest occupied energy level

  • is_spin (bool) – indicates if projected density is spin resolved

  • num_interactions (int) – number of interactions in COOP

  • interactions (list[str]) – list of the interactions included in the file

Methods

_read_coopcar(file_name='COOPCAR.lobster')[source]

Read lobster COOPCAR and extract projected overlap

Parameters

file_name (str) – file location of the COOPCAR.lobster file file

Variables

_pcoop (numpy.ndarray) – numpy array that contains the energy of levels and the projected crystal orbital overlap population densities

Returns

  • emax (float) – maximum energy level

  • emin (float) – minimum energy level

  • ndos (int) – number of descritized energy levels

  • e_fermi (float) – highest occupied energy level

  • is_spin (bool) – indicates if projected density is spin resolved

  • num_interactions (int) – number of interactions in COOP

  • interactions (list[str]) – list of the interactions included in the file

apply_gaussian_filter(sigma)[source]

Applies Gaussian filter to self._pcoop

Parameters

sigma (float) – standard deviation of Gaussian Kernel

Variables

_pcoop_original (numpy.ndarray) – self._pcoop array without filter

Notes

Understand carefully what this does before using it. It applies a Gaussian filter to the average and integrated PCOOP and the PCOOP such that the average and integrated PCOOP may become meaningless

get_average_int_pcoop(sum_spin=True)[source]

obtain average integrated overlap population for all atom pairs

Parameters

sum_spin (bool) – indicates whether data of different spins should be summed

Returns

average_int_pcoop – 1-D or 2-D array of average integrated pcoop for all interactions

Return type

numpy.ndarray

get_average_pcoop(sum_spin=True)[source]

obtain average overlap population for all atom pairs

Parameters

sum_spin (bool) – indicates whether data of different spins should be summed

Returns

average_pcoop – 1-D or 2-D array of average pcoop for all interactions

Return type

numpy.ndarray

get_bonding_states(orbital_indices, dos_energies, interactions=[], set_antibonding_zero=False, emax=inf)[source]

method for obtaining bonding fraction of dos or dos-like array

Parameters
  • orbital_indices (list[list]) – list of energy indices for each molecular orbital

  • dos_energies (numpy.ndarray) – energies at which dos is calculated

  • interactions (list) – indices of interactions for which to find the integrated pcoop

  • set_antibonding_zero (bool) – if true, set antibonding populations to zero to look at total instead of net bonding characteristics

  • emax (float) – maximum energy level

Returns

dos_bonding – 1-D array of the relative bonding for either dos-like array

Return type

numpy.ndarray

get_energies()[source]

method for obtaining energy levels

Returns

energies – 1-D array of energies

Return type

numpy.ndarray

get_integrated_pcoop(interactions=[], sum_pcoop=False, sum_spin=True, set_antibonding_zero=False)[source]

obtain integrated projected crystal orbital overlap populations

Parameters
  • interactions (list) – indices of interactions for which to find the integrated pcoop

  • sum_pcoop (bool) – indicates whether all pcoop should be summed

  • sum_spin (bool) – indicates whether data of different spins should be summed

  • set_antibonding_zero (bool) – if true, set antibonding populations to zero to look at total instead of net bonding characteristics

Returns

integrated_pcoop – 1-D, 2-D, or 3-D array of integrated pcoop for all interactions

Return type

numpy.ndarray

get_pcoop(interactions=[], sum_pcoop=False, sum_spin=True, set_antibonding_zero=False)[source]

method for obtaining projected crystal orbital overlap populations

Parameters
  • interactions (list) – indices of interactions for which to find the integrated pcoop

  • sum_pcoop (bool) – indicates whether all pcoop should be summed

  • sum_spin (bool) – indicates whether data of different spins should be summed

  • set_antibonding_zero (bool) – if true, set antibonding populations to zero to look at total instead of net bonding characteristics

Returns

pcoop – 1-D or 2-D array of pcoop for all interactions

Return type

numpy.ndarray