lumi Package¶
lumi
Package¶
deltaSCF
Module¶
- class abipy.lumi.deltaSCF.DeltaSCF(structuregs, structureex, forces_gs, forces_ex, ag_energy, ag_star_energy, ae_star_energy, ae_energy, meta=None)[source]¶
Bases:
object
Object to post-process the results from a LumiWork, following a one-effective phonon mode model (1D-CCM). For equations, notations and formalism, please refer to :
- classmethod from_json_file(json_path)[source]¶
Create the object from a json file containing the path to netcdf files, produced at the end of a LumiWork
- classmethod from_four_points_file(filepaths)[source]¶
Create the object from a list of netcdf files in the order (Ag,Agstar,Aestar,Ae). Ag: Ground state at relaxed ground state atomic positions. Agstar: Excited state at relaxed ground state atomic positions. Aestar: Excited state at relaxed excited state atomic positions. Ae: Ground state at relaxed excited state atomic positions.
- Parameters:
filepaths – list of netcdf files in the order [Ag,Agstar,Aestar,Ae]
- Returns:
A DeltaSCF object
- classmethod from_relax_file(filepaths)[source]¶
Create the object from the two relaxation files (relax_gs, relax_ex). Give only acccess to structural relaxation induced by the transition and to E_zpl
- diff_pos_mass_weighted()[source]¶
Difference between gs and ex structures in Angström, weighted by the squared atomic masses (n_atoms,3) shape
- defect_index(defect_symbol)[source]¶
Defect index in the structure from its symbol, ex defect_index(“Eu”).
- get_dataframe_atoms(defect_symbol)[source]¶
Panda dataframe with relevant properties per atom. Units : [ (mass,amu), (deltaR,Angstrom), (DeltaQ^2,amu.Angstrom^2), (DeltaF,eV/Angstrom) ]
- delta_q(unit='atomic')[source]¶
Total Delta_Q
- Parameters:
unit – amu^1/2.Angstrom if unit = ‘atomic’, kg^1/2.m if ‘SI’
- E_FC_ex(unit='eV')[source]¶
Franck condon energy in excited state (eV) = Relaxation energy between Ag* and Ae* states
- E_FC_gs(unit='eV')[source]¶
Franck condon energy in ground state (eV) = Relaxation energy between Ae and Ag states
- FWHM_1D(T=0)[source]¶
Full width at half-maximum following a semi-classical approx in the 1D-CCM (eq.20-21 of https://doi.org/10.1103/PhysRevB.96.125132)
- Parameters:
T – Temperature
- FC_factor_approx(n)[source]¶
FC factor between initial vib. state m=0 to final vib. state n approx of same eff frequency in gs and ex and T = 0K. See eq. (9) of https://doi.org/10.1002/adom.202100649
- lineshape_1D_zero_temp(energy_range=[0.5, 5], max_m=25, phonon_width=0.01, with_omega_cube=True, normalized='Area')[source]¶
Compute the emission lineshape following the effective phonon 1D-CCM at T=0K. See eq. (9) of https://doi.org/10.1002/adom.202100649.
- Parameters:
energy_range – Energy range at which the intensities are computed, ex : [0.5,5]
max_m – Maximal vibrational state m considered
phonon_width – fwhm of each phonon peak, in eV
with_omega_cube – Considered or not the omega^3 dependence of the intensity
normalized – Normalisation procedure. ‘Area’ if Area under the curve = 1. ‘Sum’ if maximum of the curve = 1.
- Returns:
E_x = Energies at which the intensities are computed I = Intensities
- plot_lineshape_1D_zero_temp(energy_range=[0.5, 5], max_m=25, phonon_width=0.01, with_omega_cube='True', normalized='Area', ax=None, **kwargs)[source]¶
Plot the the emission lineshape following the effective phonon 1D-CCM at T=0K.
- Parameters:
ax –
matplotlib.axes.Axes
or None if a new figure should be created.energy_range – Energy range at which the intensities are computed, ex : [0.5,5]
max_m – Maximal vibrational state m considered
phonon_width – fwhm of each phonon peak, in eV
with_omega_cube – Considered or not the omega^3 dependence of the intensity
normlized – Normalisation procedure. ‘Area’ if Area under the curve = 1. ‘Sum’ if maximum of the curve = 1.
Returns:
matplotlib.figure.Figure
Keyword arguments controlling the display of the figure:
kwargs
Meaning
title
Title of the plot (Default: None).
show
True to show the figure (default: True).
savefig
“abc.png” or “abc.eps” to save the figure to a file.
size_kwargs
Dictionary with options passed to fig.set_size_inches e.g. size_kwargs=dict(w=3, h=4)
tight_layout
True to call fig.tight_layout (default: False)
ax_grid
True (False) to add (remove) grid from all axes in fig. Default: None i.e. fig is left unchanged.
ax_annotate
Add labels to subplots e.g. (a), (b). Default: False
fig_close
Close figure. Default: False.
plotly
Try to convert mpl figure to plotly.
- get_dataframe(label=None)[source]¶
Panda dataframe with the main results of a LumiWork : transition energies, delta Q, Huang Rhys factor,… Units used are Angstrom, eV, amu. DeltaSCF object should be instantiated with the four points files, not with relax files only.
- draw_displacements_vesta(in_path, mass_weighted=False, scale_vector=20, width_vector=0.3, color_vector=[255, 0, 0], centered=True, factor_keep_vectors=0.1, out_path='VESTA_FILES', out_filename='gs_ex_relaxation')[source]¶
Draw the ground state to excited state atomic relaxation on a vesta structure.
- Parameters:
in_path – path where the initial .vesta structure in stored, should correspond to the ground state relaxed structure.
mass_weighted – If True, weight the displacements by the atomic masses. Draw the Delta Q in that case.
scale_vector – scaling factor of the vector modulus
width_vector – vector width
color_vector – color in rgb format
centered – center the vesta structure around [0,0,0]
factor_keep_vectors – draw only the eigenvectors with magnitude > factor_keep_vectors * max(magnitude)
out_path – path where .vesta files with vector are stored
- displacements_visu(a_g=10, **kwargs)[source]¶
Make a 3d visualisation of the displacements induced by the electronic transition = Difference between ground state and excited state atomic positions. The colors of the atoms are based on Delta_Q_^2 per atom.
- Parameters:
a_g – coefficient that multiplies the displacement magnitudes
Returns:
matplotlib.figure.Figure
Keyword arguments controlling the display of the figure:
kwargs
Meaning
title
Title of the plot (Default: None).
show
True to show the figure (default: True).
savefig
“abc.png” or “abc.eps” to save the figure to a file.
size_kwargs
Dictionary with options passed to fig.set_size_inches e.g. size_kwargs=dict(w=3, h=4)
tight_layout
True to call fig.tight_layout (default: False)
ax_grid
True (False) to add (remove) grid from all axes in fig. Default: None i.e. fig is left unchanged.
ax_annotate
Add labels to subplots e.g. (a), (b). Default: False
fig_close
Close figure. Default: False.
plotly
Try to convert mpl figure to plotly.
- plot_delta_R_distance(defect_symbol, colors=['k', 'r', 'g', 'b', 'c', 'm'], ax=None, **kwargs)[source]¶
Plot DeltaR vs distance from defect for each atom, colored by species.
- Parameters:
ax –
matplotlib.axes.Axes
or None if a new figure should be created.defect_symbol – defect_symbol, defect location will be the reference
colors – list of colors for the species
Returns:
matplotlib.figure.Figure
Keyword arguments controlling the display of the figure:
kwargs
Meaning
title
Title of the plot (Default: None).
show
True to show the figure (default: True).
savefig
“abc.png” or “abc.eps” to save the figure to a file.
size_kwargs
Dictionary with options passed to fig.set_size_inches e.g. size_kwargs=dict(w=3, h=4)
tight_layout
True to call fig.tight_layout (default: False)
ax_grid
True (False) to add (remove) grid from all axes in fig. Default: None i.e. fig is left unchanged.
ax_annotate
Add labels to subplots e.g. (a), (b). Default: False
fig_close
Close figure. Default: False.
plotly
Try to convert mpl figure to plotly.
- plot_delta_F_distance(defect_symbol, colors=['k', 'r', 'g', 'b', 'c', 'm'], ax=None, **kwargs)[source]¶
Plot DeltaF vs distance from defect for each atom, colored by species.
- Parameters:
ax –
matplotlib.axes.Axes
or None if a new figure should be created.defect_symbol – defect_symbol, defect location will be the reference
colors – list of colors for the species
Returns:
matplotlib.figure.Figure
Keyword arguments controlling the display of the figure:
kwargs
Meaning
title
Title of the plot (Default: None).
show
True to show the figure (default: True).
savefig
“abc.png” or “abc.eps” to save the figure to a file.
size_kwargs
Dictionary with options passed to fig.set_size_inches e.g. size_kwargs=dict(w=3, h=4)
tight_layout
True to call fig.tight_layout (default: False)
ax_grid
True (False) to add (remove) grid from all axes in fig. Default: None i.e. fig is left unchanged.
ax_annotate
Add labels to subplots e.g. (a), (b). Default: False
fig_close
Close figure. Default: False.
plotly
Try to convert mpl figure to plotly.
- plot_four_BandStructures(nscf_files, ax_mat=None, ylims=(-5, 5), **kwargs)[source]¶
plot the 4 band structures. nscf_files is the list of Ag, Agstar, Aestar, Ae nscf gsr file paths.
Keyword arguments controlling the display of the figure:
kwargs
Meaning
title
Title of the plot (Default: None).
show
True to show the figure (default: True).
savefig
“abc.png” or “abc.eps” to save the figure to a file.
size_kwargs
Dictionary with options passed to fig.set_size_inches e.g. size_kwargs=dict(w=3, h=4)
tight_layout
True to call fig.tight_layout (default: False)
ax_grid
True (False) to add (remove) grid from all axes in fig. Default: None i.e. fig is left unchanged.
ax_annotate
Add labels to subplots e.g. (a), (b). Default: False
fig_close
Close figure. Default: False.
plotly
Try to convert mpl figure to plotly.
- draw_displaced_parabolas(ax=None, scale_eff_freq=4, font_size=8, **kwargs)[source]¶
Draw the four points diagram with relevant transition energies.
Args: ax:
matplotlib.axes.Axes
or None if a new figure should be created. scale_eff_freq: scaling factor to adjust the parabolas curvatures. font_size: font size for the annotationsReturns:
matplotlib.figure.Figure
Keyword arguments controlling the display of the figure:
kwargs
Meaning
title
Title of the plot (Default: None).
show
True to show the figure (default: True).
savefig
“abc.png” or “abc.eps” to save the figure to a file.
size_kwargs
Dictionary with options passed to fig.set_size_inches e.g. size_kwargs=dict(w=3, h=4)
tight_layout
True to call fig.tight_layout (default: False)
ax_grid
True (False) to add (remove) grid from all axes in fig. Default: None i.e. fig is left unchanged.
ax_annotate
Add labels to subplots e.g. (a), (b). Default: False
fig_close
Close figure. Default: False.
plotly
Try to convert mpl figure to plotly.
lineshape
Module¶
- class abipy.lumi.lineshape.Lineshape(E_zpl, ph_eigvec, ph_eigfreq, structure, forces, displacements, use_forces)[source]¶
Bases:
object
Object representing a luminescent lineshape, following a multi-phonon mode model (multiD-CCM). For 1D-CCM, use plot_lineshape_1D_zero_temp() function of the abipy/lumi/deltaSCF module.
- For equations, notations and formalism, please refer to :
https://doi.org/10.1103/PhysRevB.96.125132 https://doi.org/10.1002/adom.202100649 https://pubs.acs.org/doi/full/10.1021/acs.chemmater.3c00537
In the 1D-CCM, the vibronic peaks are the ones from a fictious phonon mode that connects the atomic relaxation between the ground state and excited state. Within this model, the global shape (fwhm) is well represented if the total Huang-Rhys factor is large enough (gaussian shaped spectrum). However, if vibronic peaks are present in the experimental spectrum, this model is not able to correctly reproduce them as it assumes a fictious phonon mode.
In the multiD-CCM, the atomic relaxation is projected along the phonon eigenvectors of the system, allowing a phonon-projected decomposition of the relaxation. Better agreement with experimental vibronic peaks is expected.
- classmethod from_phonopy_phonons(E_zpl, phonopy_ph, dSCF_structure, use_forces=True, dSCF_displacements=None, dSCF_forces=None, coords_defect_dSCF=None, tol=0.3)[source]¶
Different levels of approximations for the phonons and force/displacements: See discussion in the supplementary informations of https://pubs.acs.org/doi/full/10.1021/acs.chemmater.3c00537, section (1).
size_supercell deltaSCF = size_supercell phonons (phonons of the bulk structure or phonons of defect structure). Use of the forces or the displacemements is allowed.
size_supercell dSCF < size_bulk_supercell phonons (bulk) Use of the forces only.
size_supercell dSCF < size__defect_supercell phonons (embedding) Use of the forces only
The code first extracts the eigenmodes of the phonopy object. Then, it tries performs a structure matching between the phonopy structure and the dSCF structure (critical part) in order to put the displacements/forces on the right atoms.
- Parameters:
E_zpl – Zero-phonon line energy in eV
phonopy_ph – Phonopy object containing eigenfrequencies and eigenvectors
dSCF_structure – Delta_SCF structure
dSCF_displacements – Dispalcements Delta R induced by the electronic transition
dSCF_forces – Dispalcements Delta F induced by the electronic transition
coords_defect_dSCF – Main coordinates of the defect in defect structure, if defect complex, can be set to the center of mass of the complex
tol – tolerance in Angstrom applied for the matching between the dSCF structure and phonon structure
Returns: A lineshape object
- p_nu()[source]¶
Contribution of each mode, eq. (11) of https://doi.org/10.1002/adom.202100649
- eff_freq_multiD()[source]¶
Effective coupling frequency, eq. (13) of https://doi.org/10.1002/adom.202100649
- S_hbarOmega(broadening)[source]¶
Return (Energy,spectral decomposition S(hw))
- Parameters:
broadening – fwhm of the gaussian broadening in meV
- Bose_einstein(T, freq)[source]¶
Bose Einstein average occupation number
- Parameters:
T – Temperature in K
freq – frequency in eV
- G_t(T, S_nu, omega_nu)[source]¶
Generation function. Eq. (3) of https://pubs.acs.org/doi/full/10.1021/acs.chemmater.3c00537
- Parameters:
T – Temperature in K
S_nu – Parial Huang-Rhys factors
omega_nu – phonon frequencies
- A_hw(T, lamb=5, w=1, model='multi-D')[source]¶
Lineshape function Eq. (2) of https://pubs.acs.org/doi/full/10.1021/acs.chemmater.3c00537 Returns (Energy in eV, Lineshape function )
- Parameters:
T – Temperature in K
lamb – Lorentzian broadening applied to the vibronic peaks, in meV
w – Gaussian broadening applied to the vibronic peaks, in meV
model – ‘multi-D’ for full phonon decomposition, ‘one-D’ for 1D-CCM PL spectrum.
- L_hw(T=0, lamb=5, w=1, model='multi-D')[source]¶
Normalized Luminescence intensity (area under the curve = 1) Eq. (1) of https://pubs.acs.org/doi/full/10.1021/acs.chemmater.3c00537 Returns (Energy in eV, Luminescence intensity)
- Parameters:
T – Temperature in K
lamb – Lorentzian broadening applied to the vibronic peaks, in meV
w – Gaussian broadening applied to the vibronic peaks, in meV
model – ‘multi-D’ for full phonon decomposition, ‘one-D’ for 1D-CCM PL spectrum.
- plot_spectral_function(broadening=1, ax=None, with_S_nu=False, **kwargs)[source]¶
Plot the Huang-Rhys spectral function S_hbarOmega
- Parameters:
broadening – fwhm of the gaussian broadening in meV
with_S_nu – True to add stem lines associated to the individuals partial Huang-Rhys factors
Keyword arguments controlling the display of the figure:
kwargs
Meaning
title
Title of the plot (Default: None).
show
True to show the figure (default: True).
savefig
“abc.png” or “abc.eps” to save the figure to a file.
size_kwargs
Dictionary with options passed to fig.set_size_inches e.g. size_kwargs=dict(w=3, h=4)
tight_layout
True to call fig.tight_layout (default: False)
ax_grid
True (False) to add (remove) grid from all axes in fig. Default: None i.e. fig is left unchanged.
ax_annotate
Add labels to subplots e.g. (a), (b). Default: False
fig_close
Close figure. Default: False.
plotly
Try to convert mpl figure to plotly.
- plot_emission_spectrum(unit='eV', T=0, lamb=3, w=3, ax=None, **kwargs)[source]¶
Plot the Luminescence intensity
- Parameters:
unit – ‘eV’, ‘cm-1’, or ‘nm’
T – Temperature in K
lamb – Lorentzian broadening applied to the vibronic peaks, in meV
w – Gaussian broadening applied to the vibronic peaks, in meV
Keyword arguments controlling the display of the figure:
kwargs
Meaning
title
Title of the plot (Default: None).
show
True to show the figure (default: True).
savefig
“abc.png” or “abc.eps” to save the figure to a file.
size_kwargs
Dictionary with options passed to fig.set_size_inches e.g. size_kwargs=dict(w=3, h=4)
tight_layout
True to call fig.tight_layout (default: False)
ax_grid
True (False) to add (remove) grid from all axes in fig. Default: None i.e. fig is left unchanged.
ax_annotate
Add labels to subplots e.g. (a), (b). Default: False
fig_close
Close figure. Default: False.
plotly
Try to convert mpl figure to plotly.
- abipy.lumi.lineshape.get_forces_on_phonon_supercell(dSCF_supercell, phonon_supercell, forces_dSCF, tol)[source]¶