pytaser package#

pytaser.generator module#

This module contains the TASGenerator class, which is used to generate TAS spectra.

class pytaser.generator.TASGenerator(bs, kpoint_weights, dos, dfc=None)[source]#

Bases: object

Class to generate a TAS spectrum (decomposed and cumulative) from a bandstructure and dos object.

Parameters:
  • bs – Pymatgen-based bandstructure object

  • kpoint_weights – kpoint weights either found by the function or inputted.

  • dos – Pymatgen-based dos object

  • dfc – Pymatgen-based DielectricFunctionCalculator object (for computing oscillator strengths).

bs#

Pymatgen bandstructure object

kpoint_weights#

k-point weights (degeneracies).

dos#

Pymatgen-based dos object

dfc#

Pymatgen-based DielectricFunctionCalculator object (for computing oscillator strengths)

bg_centre#

Energy (eV) of the bandgap centre.

vb#

Spin dict detailing the valence band maxima.

cb#

Spin dict detailing the conduction band minima

band_occupancies(temp, conc, dark=True)[source]#

Gives band occupancies.

Parameters:
  • temp – Temperature of material we wish to investigate (affects the FD distribution)

  • conc – Carrier concentration of holes and electrons (both are the same). Inversely proportional to pump-probe time delay.

  • dark – Bool; dark = True indicates pump is on.

Returns:

occs} for all bands across all k-points.

Return type:

A dictionary of {Spin

classmethod from_mpid(mpid, bg=None, api_key=None, mpr=None)[source]#

Import the desired bandstructure and dos objects from the legacy Materials Project database.

Parameters:
  • mpid – The Materials Project ID of the desired material.

  • bg – The experimental bandgap (eV) of the material. If None, the band gap of the MP calculation will be used.

  • api_key – The user’s Materials Project API key.

  • mpr – An MPRester object if already generated by user.

Returns:

A TASGenerator object.

classmethod from_vasp_outputs(vasprun_file, waveder_file=None, bg=None)[source]#

Create a TASGenerator object from VASP output files.

Parameters:
  • vasprun_file – Path to vasprun.xml file (to generate bandstructure object).

  • waveder_file – Path to WAVEDER file (to generate dielectric function calculator object,

  • moments). (to compute orbital derivatives and transition dipole) –

  • bg – Experimental bandgap of the material, used to rigidly shift the DFT results to match

  • None (this value. If) –

Returns:

A TASGenerator object.

generate_tas(temp, conc, energy_min=0, energy_max=5, gaussian_width=0.1, cshift=None, step=0.01, light_occs=None, dark_occs=None, processes=None)[source]#

Generates TAS spectra based on inputted occupancies, and a specified energy mesh. If the TASGenerator has not been generated from VASP outputs (and thus does not have a dfc attribute), then the output TAS is generated using the change in joint density of states (JDOS) under illumination, with no consideration of oscillator strengths. Otherwise, the output TAS is generated considering all contributions to the predicted TAS spectrum.

Parameters:
  • temp – Temperature (K) of material we wish to investigate (affects the FD distribution)

  • conc – Carrier concentration (cm^-3) of holes and electrons (both are equivalent). Inversely proportional to pump-probe time delay.

  • energy_min – Minimum band transition energy to consider for energy mesh (eV)

  • energy_max – Maximum band transition energy to consider for energy mesh (eV)

  • gaussian_width – Width of gaussian curve

  • cshift – Complex shift in the Kramers-Kronig transformation of the dielectric function (see https://www.vasp.at/wiki/index.php/CSHIFT). If not set, uses the value of CSHIFT from the underlying VASP WAVEDER calculation. (only relevant if the TASGenerator has been generated from VASP outputs)

  • step – Interval between energy points in the energy mesh.

  • light_occs – Optional input parameter for occupancies of material under light, otherwise automatically calculated based on input temperature (temp) and carrier concentration (conc) [dict]

  • dark_occs – Optional input parameter for occupancies of material in dark, otherwise automatically calculated based on input temperature (temp) and carrier concentration (conc) [dict]

  • processes – Number of processes to use for multiprocessing. If not set, defaults to one less than the number of CPUs available.

Returns:

TAS class containing the following inputs;
  • tas_total: overall deltaT TAS spectrum for a material under the

    specified conditions.

  • jdos_diff_if: JDOS difference (from dark to light) across the energy mesh for a

    specific band transition i (initial) -> f (final) [dict]

  • jdos_light_total: overall JDOS (pump-on) for a material under the

    specified conditions

  • jdos_light_if: JDOS (pump-on) across the energy mesh for a specific band

    transition i (initial) -> f (final) [dict]

  • jdos_dark_total: overall JDOS (pump-off) for a material under the

    specified conditions

  • jdos_dark_if: JDOS (pump-off) across the energy mesh for a specific band

    transition i (initial) -> f (final) [dict]

  • energy_mesh_ev: Energy mesh of spectra in eV, with an interval of ‘step’.

  • bandgap: Bandgap of the system, in eV, rounded to 2 decimal points

  • temp: Temperature of the system, in K

  • conc: Carrier concentration of the system, in cm^-3

  • alpha_dark: Absorption coefficient of the material in the dark, in cm^-1 (only

    calculated if the TASGenerator has been generated from VASP outputs)

  • alpha_light_dict: Dictionary of band-to-band absorption, stimulated emission

    and summed contributions to the total overall absorption coefficient under illumination, in cm^-1 (only calculated if the TASGenerator has been generated from VASP outputs)

  • weighted_jdos_diff_if: JDOS difference (from dark to light) across the energy

    mesh for a specific band transition i (initial) -> f (final), weighted by the oscillator strength of the transition [dict]

  • weighted_jdos_light_if: JDOS (pump-on) across the energy mesh for a specific band

    transition i (initial) -> f (final), weighted by the oscillator strength of the transition [dict]

  • weighted_jdos_dark_if: JDOS (pump-off) across the energy mesh for a specific band

    transition i (initial) -> f (final), weighted by the oscillator strength of the transition [dict]

pytaser.generator.gaussian(x, width, center=0.0, height=None)[source]#

Returns Gaussian curve(s) centred at point(s) x, where x is array-like.

Parameters:
  • x – Input array.

  • width – Width of the gaussian.

  • center – Center of the gaussian.

  • height – height of the gaussian. If height is None, a normalized gaussian is returned.

pytaser.generator.get_cbm_vbm_index(bs)[source]#
Parameters:

bs – bandstructure object.

Returns:

Valence and Conduction band as an index number for different spins.

pytaser.generator.get_nonzero_band_transitions(occs, eigs_shifted, rspin, spin, sigma, nedos, deltae, ismear, min_band, max_band, nk)[source]#

Helper function to filter band transitions before (multi)processing.

pytaser.generator.jdos(bs, f, i, occs, energies, kweights, gaussian_width, spin=Spin.up)[source]#

Obtains the cumulative JDOS value for a specific i->f transition, with consideration of partial occupancy and spin polarisation.

Parameters:
  • bs – bandstructure object

  • f – final band

  • i – initial band

  • occs – occupancies over all bands.

  • energies – energy mesh (eV)

  • kweights – k-point weights

  • gaussian_width – width of gaussian plot.

  • spin – Which spin channel to include.

Returns:

Cumulative JDOS value for a specific i->f transition, with consideration of partial occupancy and spin polarisation.

pytaser.generator.occ_dependent_alpha(dfc, occs, spin=Spin.up, sigma=None, cshift=None, processes=None, energy_max=6)[source]#

Calculate the expected optical absorption given the groundstate orbital derivatives and eigenvalues (via dfc) and specified band occupancies. Templated from pymatgen.io.vasp.optics.epsilon_imag().

Parameters:
  • dfc – DielectricFunctionCalculator object

  • occs – Array of band occupancies with shape (nspin, nbands, nkpoints)

  • spin – Which spin channel to include.

  • sigma – Smearing width (in eV) for broadening of the dielectric function (see https://www.vasp.at/wiki/index.php/SIGMA). If not set, uses the value of SIGMA from the underlying VASP WAVEDER calculation.

  • cshift – Complex shift in the Kramers-Kronig transformation of the dielectric function (see https://www.vasp.at/wiki/index.php/CSHIFT). If not set, uses the value of CSHIFT from the underlying VASP WAVEDER calculation.

  • processes – Number of processes to use for multiprocessing. If not set, defaults to one less than the number of CPUs available.

  • energy_max – Maximum band transition energy to consider (in eV). A minimum range of -6 eV to +6 eV is considered regardless of energy_max, to ensure a reasonable dielectric function output.

Returns:

(alpha_dict, tdm_array) where alpha_dict is a dictionary of band-to-band absorption, stimulated emission and summed contributions to the total overall absorption coefficient under illumination in cm^-1, and tdm_array is an array of shape (nbands, nbands, nkpoints) with the transition dipole matrix elements for each band-to-band transition.

pytaser.generator.set_bandgap(bandstructure, dos, bandgap)[source]#

Shifts all bands of a material to correct the DFT-underestimated bandgap according to the input experimental bandgap.

Parameters:
  • bandstructure – PMG bandstructure object

  • dos – PMG dos object

  • bandgap – experimental bandgap from literature (float)

Returns:

new bandstructure and dos objects for the same material, but with corrected bandgap.

pytaser.das_generator module#

Created on Thu Aug 3 16:42:36 2023.

@author: lucasverga

class pytaser.das_generator.DASGenerator(new_system, reference_system)[source]#

Bases: object

Class to generate a DAS spectrum (decomposed and cumulative) from a bandstructure and dos object.

Parameters:
  • new_system – Internal_Abs object from internal_abs_generator for the new system

  • reference_system – Internal_Abs object from internal_abs_generator for the reference system

new_system#

Internal_Abs object from internal_abs_generator for the new system

reference_system#

Internal_Abs object from internal_abs_generator for the reference system

Parameters:
  • new_system – Internal_Abs object from internal_abs_generator for the new system

  • reference_system – Internal_Abs object from internal_abs_generator for the reference system.

classmethod from_mpid(mpid, mpid_ref, bg=None, bg_ref=None, api_key=None, mpr=None, mpr_ref=None)[source]#

Import the desired bandstructure and dos objects from the legacy Materials Project database.

Parameters:
  • mpid – The Materials Project ID of the new system.

  • mpid_ref – The Materials Project ID of the reference system.

  • bg – The experimental bandgap (eV) of the new system. If None, the band gap of the MP calculation will be used.

  • bg_ref – The experimental bandgap (eV) of the reference system. If None, the band gap of the MP calculation will be used.

  • api_key – The user’s Materials Project API key.

  • mpr – An MPRester object for the new system if already generated by user.

  • mpr_ref – An MPRester object for the reference system if already generated by user.

Returns:

A DASGenerator object containing the Internal_Abs object for the new system and reference system.

classmethod from_vasp_outputs(vasprun_file_new_system, vasprun_file_ref, waveder_file_new_system=None, waveder_file_ref=None)[source]#

Create a DASGenerator object from VASP output files.

The user should provide the vasprun files for the new system and the reference system, followed by the waveder files for the new system and the reference system.

Parameters:
  • vasprun_file_new_system – The vasprun.xml file for the new system.

  • vasprun_file_ref – The vasprun.xml file for the reference system.

  • waveder_file_new_system – The WAVEDER file for the new system.

  • waveder_file_ref – The WAVEDER file for the reference system.

Returns:

A DASGenerator object containing the Internal_Abs object for the new system and reference system.

generate_das(temp, energy_min=0, energy_max=5, gaussian_width=0.1, cshift=None, step=0.01, new_sys_occs=None, ref_occs=None, processes=None)[source]#

Generates DAS spectra (new system - reference system) based on inputted occupancies, and a specified energy mesh. If the DASGenerator has not been generated from VASP outputs (and thus does not have a dfc attribute), then the output DAS is generated using the change in joint density of states (JDOS) from both systems, with no consideration of oscillator strengths. Otherwise, the output DAS is generated considering all contributions to the predicted DAS spectrum.

Parameters:
  • temp – Temperature (K) of material we wish to investigate (affects the FD distribution)

  • energy_min – Minimum band transition energy to consider for energy mesh (eV)

  • energy_max – Maximum band transition energy to consider for energy mesh (eV)

  • gaussian_width – Width of gaussian curve

  • cshift – Complex shift in the Kramers-Kronig transformation of the dielectric function (see https://www.vasp.at/wiki/index.php/CSHIFT). If not set, uses the value of CSHIFT from the underlying VASP WAVEDER calculation. (only relevant if the DASGenerator has been generated from VASP outputs)

  • step – Interval between energy points in the energy mesh.

  • new_sys_occs – Optional input parameter for occupancies of the new system, otherwise automatically calculated based on input temperature (temp)

  • ref_occs – Optional input parameter for occupancies of the reference system, otherwise automatically calculated based on input temperature (temp)

  • processes – Number of processes to use for multiprocessing. If not set, defaults to one less than the number of CPUs available.

Returns:

DAS class containing the following inputs;
  • das_total: overall deltaT DAS spectrum for new system - reference system.

  • jdos_new_sys_total: overall JDOS for the new system.

  • jdos_new_sys_if: JDOS for the new system across the energy mesh for a specific band

    transition i (initial) -> f (final) [dict]

  • jdos_ref_total: overall JDOS for the reference system.

  • jdos_ref_if: JDOS for the reference system across the energy mesh for a specific band

    transition i (initial) -> f (final) [dict]

  • energy_mesh_ev: Energy mesh of spectra in eV, with an interval of ‘step’.

  • bandgap: Bandgap of the system, in eV, rounded to 2 decimal points

  • temp: Temperature of the system, in K

  • alpha_ref: Absorption coefficient of the reference system, in cm^-1 (only

    calculated if the DASGenerator has been generated from VASP outputs)

  • alpha_new_sys: Absorption coefficient of the new system, in cm^-1 (only

    calculated if the DASGenerator has been generated from VASP outputs

  • weighted_jdos_diff_if: JDOS difference (from reference to new system) across the energy

    mesh for a specific band transition i (initial) -> f (final), weighted by the oscillator strength of the transition [dict]

  • weighted_jdos_new_sys_if: JDOS of new system across the energy mesh for a specific band

    transition i (initial) -> f (final), weighted by the oscillator strength of the transition [dict]

pytaser.kpoints module#

This module generates kpoint-weights for uniform-mesh non-magnetic materials.

This is vital when using the Materials Project database to generate spectra in PyTASER.

pytaser.kpoints.expand_kpoints(structure, kpoints, symprec=0.01, return_mapping=False, time_reversal=True)[source]#

Function to expand the kpoints.

Parameters:
  • structure – PMG structure object

  • kpoints – uniform mesh kpoints array.

  • symprec – Symmetry precision in Angstrom.(Lower value is more precise, but computationally more expensive)

  • return_mapping – bool

  • time_reversal – Time reversal operator, bool (True for +1, False for -1)

Returns:

full_kpoints, rotations, translations, is_tr, op_mapping, kp_mapping

pytaser.kpoints.get_kpoint_weights(bandstructure, time_reversal=True, symprec=0.1)[source]#

Function to calculate the kpoint_weights for non-magnetic materials (non- metals).

Parameters:
  • bandstructure – PMG bandstructure object

  • time_reversal – Time reversal operator, bool (True for +1, False for -1)

  • symprec – Symmetry precision in Angstrom.(Lower value is more precise, but computationally more expensive)

Returns:

k-point_weights

pytaser.kpoints.get_kpoints_from_bandstructure(bandstructure, cartesian=False)[source]#

Function to pull the kpoint from the bandstructure.

Parameters:
  • bandstructure – PMG bandstructure object

  • cartesian – bool, indicate if cartesian or fractional coordinates used.

Returns:

k-points

pytaser.kpoints.get_mesh_from_kpoint_diff(kpoints, ktol=1e-05)[source]#

Function to get the uniform mesh from kpoint differences.

Parameters:
  • kpoints – uniform mesh kpoints array.

  • ktol – threshold for filtering changes in kpoint-mesh points.

Returns:

np.array([na, nb, nc]), is_shifted

pytaser.kpoints.get_reciprocal_point_group_operations(structure: Structure, symprec: float = 0.01, time_reversal: bool = True)[source]#

Function to get the reciprocal point group operations.

Parameters:
  • structure – PMG structure object

  • symprec – Symmetry precision in Angstrom.(Lower value is more precise, but computationally more expensive)

  • time_reversal – Time reversal operator, bool (True for +1, False for -1)

Returns:

rotations[sort_idx], translations[sort_idx], is_tr[sort_idx]

pytaser.plotter module#

TAS Plotting Tools.

class pytaser.plotter.TASPlotter(container, material_name=None, system_name=None, reference_name=None)[source]#

Bases: object

Class to generate a matplotlib plot of the TAS, DAS, or JDOS spectra, with a specific energy mesh, material, and conditions.

Parameters:
  • container – TAS container class as generated by TASGenerator().generate_tas() or DASGenerator().generate_das()

  • material_name – Name of material being investigated (string) [optional, for labelling]

  • system_name – Name of the new system used for DAS (string) [optional, for labelling]

  • reference_name – Name of the reference system used for DAS (string) [optional, for labelling].

get_plot(relevant_transitions='auto', xaxis='wavelength', transition_cutoff=0.03, xmin=None, xmax=None, ymin=None, ymax=None, yaxis='tas', invert_axis=False, **kwargs)[source]#

Plots TAS spectra using the data generated in the TAS container class as generated by TASGenerator().generate_tas() or in the DAS container class generated by DASGenerator().generate_das(). If the TASGenerator was not generated from VASP outputs, then the output TAS is generated using the change in joint density of states (JDOS) under.

illumination, with no consideration of oscillator strengths (“Total TAS (ΔJDOS only)”) - this behaviour can also be explicitly chosen by setting “yaxis” to “jdos_diff”.

Otherwise, the output TAS is generated considering all contributions to the predicted TAS spectrum (“Total TAS (Δɑ)”). If only the contribution from the change in absorption is desired, with stimulated emission neglected, set “yaxis” to “tas_absorption_only”.

One can also plot the JDOS before and after illumination, by setting “yaxis” to “jdos”, or the effective absorption coefficient (a, including absorption and stimulated emission) before and after illumination, by setting “yaxis” to “alpha”.

If a DAS container is parsed one can opt to set “yaxis” to “das” or “jdos_das”. For the “yaxis=das” the difference between the optical absorption spectra of the new system and the reference system is plotted, together with the absorption spectra of both systems. This option can only be used if the DASGenerator was created using the WAVEDER VASP output. For the “yaxis=jdos_das” the difference between the JDOS of the new system and the reference system is ploteed, together with the JDOS of both systems.

The individual band-to-band transitions which contribute to the JDOS/TAS are also shown (which can be controlled with the transition_cutoff argument - set to 1 to not show any band-band transitions), and these are weighted by the oscillator strength of the transition when TASGenerator was created from VASP objects and yaxis=”tas”, “tas_absorption_only” or “alpha”.

Parameters:
  • relevant_transitions – List containing individual transitions to be displayed in the plot alongside the total plot. If material is not spin-polarised, only write the bands involved [(-1,6),(2,7),(-8,-5)…] If spin-polarised, include the type of spin involved in transition [(-1,6, “down”),(2,7, “down”),(-8,-5, “up”)…] Default is ‘auto’ mode, which shows all band transitions above the transition_cutoff. Note that band-band transitions with overlapping extrema are scaled by 95% to avoid overlapping lines.

  • xaxis – Units for the energy mesh. Either “wavelength” or “energy”.

  • transition_cutoff – The minimum height of the band transition as a percentage threshold of the largest contributing transition (Across all kpoints, within the energy range specified). Default is 0.03 (3%).

  • xmin – Minimum energy point in mesh (float64)

  • xmax – Maximum energy point in mesh (float64)

  • ymin – Minimum absorption point. Default is 1.15 * minimum point.

  • ymax – Maximum absorption point. Default is 1.15 * maximum point.

  • yaxis

    What spectral data to plot. If yaxis = “tas” (default), will plot the predicted TAS spectrum considering all contributions to the TAS (if TASGenerator was created from VASP outputs)(“Total TAS (Δɑ)”), or just using the change in the joint density of states (JDOS) before & after illumination (“Total TAS (ΔJDOS only)”) - the latter can also be explicitly selected with yaxis = “jdos_diff”.

    Alternative choices:

    • yaxis = "tas_absorption_only":

      Plot the predicted TAS spectrum considering only the change in absorption, with stimulated emission neglected.

    • yaxis = "jdos":

      Plot the joint density of states before & after illumination.

    • yaxis = "alpha":

      Plot the effective absorption coefficient (a, including absorption and stimulated emission) before and after illumination.

    • yaxis = "das":

      Plot the difference between the optical absorption spectra of the new system and the reference system, together with the absorption spectra

    • yaxis = "jdos_das":

      Plot the difference between the JDOS of the new system and the reference system, together with the JDOS of both systems.

  • invert_axis – If True, x-axis decreases from left-right. If False (default), x-axis will increase from left-right. Especially relevant when switching between wavelength and electronvolts.

  • **kwargs – Additional arguments to be passed to matplotlib.pyplot.legend(); such as ncols (number of columns in the legend), loc (location of the legend), fontsize etc. (see https://matplotlib.org/stable/api/_as_gen/matplotlib.pyplot.legend.html)

Returns:

Matplotlib pyplot of the desired spectrum, with labelled units.

pytaser.plotter.cutoff_transitions(dictionary, cutoff, ind_xmin, ind_xmax)[source]#

Output a list of transitions from a dict, with any that fall below a percentage cutoff of the maximum value transition set to None.

pytaser.plotter.ev_to_lambda(ev)[source]#

Convert photon energies from eV to a wavelength in nm.

pytaser.plotter.lambda_to_ev(lambda_float)[source]#

Convert photon energies from a wavelength in nm to eV.