hylight package

Subpackages

Submodules

hylight.constants module

A collection of useful physics constants.

hylight.guess_width module

Semi classical guess of linewidth.

class hylight.guess_width.OmegaEff(wm)

Bases: object

Mode of computation of a effective frequency.

FC_MEAN:

\[\Omega = \frac{\sum_j \omega_j d^\text{FC}_j}{\sum_j d^\text{FC}_j}\]

Should be used with WidthModel.ONED because it is associated with the idea that all the directions are softened equally in the excited state.

HR_MEAN:

\[\Omega = \frac{d^\text{FC}}{S_\text{tot}} = \frac{\sum_j \omega_j S_j}{\sum_j S_j}\]

HR_RMS:

\[\Omega = \sqrt{\frac{\sum_j \omega_j^2 S_j}{\sum_j S_j}}\]

FC_RMS:

\[\Omega = \sqrt{\frac{\sum_j \omega_j^2 d^\text{FC}_j}{\sum_j d^\text{FC}_j}}\]

Should be used with WidthModel.SINGLE_ES_FREQ because it makes sense when we get only one Omega_eff for the excited state (this single effective frequency should be computed beforehand)

GRAD:

The effective frequency of GS is computed in the direction of the gradiant of GS PES at the ES position. The effective frequency for the ES is provided by the user (and should be computed in the same direction). Should be used with WidthModel.ONED.

DISP:

The effective frequency of GS is computed in the direction of the displacement. The effective frequency for the ES is provided by the user (and should be computed in the same direction). Should be used with WidthModel.ONED.

DISP: OmegaEff = <hylight.guess_width.OmegaEff object>
FC_MEAN: OmegaEff = <hylight.guess_width.OmegaEff object>
FC_RMS: OmegaEff = <hylight.guess_width.OmegaEff object>
GRAD: OmegaEff = <hylight.guess_width.OmegaEff object>
HR_MEAN: OmegaEff = <hylight.guess_width.OmegaEff object>
HR_RMS: OmegaEff = <hylight.guess_width.OmegaEff object>
class hylight.guess_width.WidthModel(value, names=<not given>, *values, module=None, qualname=None, type=None, start=1, boundary=None)

Bases: Enum

Mode of approximation of the band width.

ONED: The witdth is computed in a 1D mode for both ES and GS.

SINGLE_ES_FREQ: We suppose all modes of the ES have the same frequency (that should be provided).

FULL_ND: (Not implemented) We know the frequency of each ES mode.

FULL_ND = 2
ONED = 0
SINGLE_ES_FREQ = 1
hylight.guess_width.duschinsky(phonons_a, phonons_b)

Dushinsky matrix from b to a \(S_{a \gets b}\).

hylight.guess_width.effective_phonon_energy(omega_eff_type, hrs, es, masses, *, delta_R=None, modes=None, mask=None)

Compute an effective phonon energy in eV following the strategy of omega_eff_type.

Parameters:
  • omega_eff_type – The mode of evaluation of the effective phonon energy.

  • hrs – The array of Huang-Rhyes factor for each mode.

  • es – The array of phonon energy in eV.

  • masses – The array of atomic masses in atomic mass unit.

  • delta_R – The displacement between GS and ES in A. It is only required if omega_eff_type is ONED_FREQ.

Returns:

The effective energy in eV.

hylight.guess_width.expected_width(phonons, delta_R, fc_shift_gs, fc_shift_es, T, mask=None, omega_eff_type=<hylight.guess_width.OmegaEff object>, width_model=WidthModel.ONED)

Compute a spectrum without free parameters.

Parameters:
  • phonons – list of modes

  • delta_R – displacement in A

  • fc_shift_gs – Ground state/absorption Franck-Condon shift in eV

  • fc_shift_es – Exceited state/emmission Franck-Condon shift in eV

  • T – temperature in K

  • resolution_e – energy resolution in eV

  • bias – ignore low energy vibrations under bias in eV

  • window_fn – windowing function in the form provided by numpy (see numpy.hamming)

  • pre_convolve – (float, optional, None) if not None, standard deviation of the pre convolution gaussian

  • shape – ZPL line shape.

  • omega_eff_type – mode of evaluation of effective frequency.

  • result_store – a dictionary to store some intermediate results.

  • ex_pes – mode of evaluation of the ES PES curvature.

Returns:

(sig(T=0), sig(T))

hylight.guess_width.gradiant_at(modes, masses, delta_R, *, mask=None)

Compute the gradient of the PES at the position of delta_R.

Parameters:
  • modes – list of modes

  • masses – list of atomic masses

  • delta_R – displacement of each atoms in A

  • mask – (optional, None) a hylight.mode.Mask used to discard irrelevant modes

Returns:

the gradient of the PES at the position of delta_R

hylight.guess_width.guess_width(phonons, delta_R, fc_shift_gs, fc_shift_es, T, resolution_e=0.001, mask=None, shape=LineShape.GAUSSIAN, omega_eff_type=<hylight.guess_width.OmegaEff object>, width_model=WidthModel.ONED, window_fn=<function hamming>, trial_line_sig=0.02)

Try to guess the width of the line from a 1D semi-classical model.

Parameters:
  • phonons – list of modes (see load_phonons()) or path to load the modes

  • delta_R – displacement in A in a numpy array (see compute_delta_R()) or tuple of two paths (pos_gs, pos_es)

  • fc_shift_gs – Ground state/absorption Franck-Condon shift in eV

  • fc_shift_es – Exceited state/emmission Franck-Condon shift in eV

  • T – temperature in K

  • resolution_e – energy resolution in eV

  • mask – a mask used in other computations to show on the plot.

  • shape – the lineshape (a LineShape instance)

  • omega_eff_type – mode of evaluation of effective frequency.

  • ex_pes – mode of evaluation of the ES PES curvature.

  • window_fn – windowing function in the form provided by numpy (see numpy.hamming)

Returns:

a guess width in eV

hylight.guess_width.integrate(x, y)

Integrate a function over x.

Parameters:
  • xnumpy.ndarray of x values

  • yy = f(x)

Returns:

\(\int f(x) dx\)

hylight.guess_width.sigma(T, S_em, e_phonon_g, e_phonon_e)

Temperature dependant standard deviation of the lineshape.

Parameters:
  • T – temperature in K

  • S_em – emmission Huang-Rhys factor

  • e_phonon_g – energy of the GS PES vibration (eV)

  • e_phonon_e – energy of the ES PES vibration (eV)

hylight.guess_width.sigma_full_nd(T, delta_R, modes_gs, modes_es, mask=None)

Compute the width of the ZPL for the ExPES.FULL_ND mode.

Parameters:
  • T – temperature in K.

  • delta_R – distorsion in A.

  • modes_gs – list of Modes of the ground state.

  • modes_es – list of Modes of the excited state.

Returns:

numpy.ndarray with only the width for the modes that are real in ground state.

hylight.guess_width.sigma_hybrid(T, S, e_phonon, e_phonon_e)

Compute the width of the ZPL for the ExPES.SINGLE_ES_FREQ mode.

hylight.guess_width.variance(x, y)

Compute the variance of a random variable of distribution y.

hylight.loader module

Wrapper for the npz loader.

hylight.mode module

Vibrational mode and related utilities.

class hylight.mode.CellMismatch(reason, details)

Bases: object

A falsy value explaining how the cell are not matching.

class hylight.mode.Mask(intervals: list[tuple[float, float]])

Bases: object

An energy based mask for the set of modes.

accept(value: float) bool

Return True if value is not under the mask.

add_interval(interval: tuple[float, float]) None

Add a new interval to the mask.

as_bool(ener: ndarray) ndarray

Convert to a boolean np.ndarray based on ener.

classmethod from_bias(bias: float) Mask

Create a mask that reject modes of energy between 0 and bias.

Parameters:

bias – minimum of accepted energy (eV)

Returns:

a fresh instance of Mask.

plot(ax, unit)

Add a graphical representation of the mask to a plot.

Parameters:
  • ax – a matplotlib Axes object.

  • unit – the unit of energy to use (ex: hylight.constant.eV_in_J if the plot uses eV)

Returns:

a function that must be called without arguments after resizing the plot.

reject(value: float) bool

Return True if value is under the mask.

class hylight.mode.Mode(lattice: ndarray, atoms: list[str], n: int, real: bool, energy: float, ref: ndarray, eigenvector: ndarray, masses: Iterable[float])

Bases: object

The representation of a vibrational mode.

It stores the eigenvector and eigendisplacement. It can be used to project other displacement on the eigenvector.

energies()

Return the energy participation of each atom to the mode.

property energy_cm1

Energy of the mode in \(cm^{-1}\).

property energy_eV

Energy of the mode in eV.

property energy_meV

Energy of the mode in meV.

property energy_si

Energy of the mode in J.

huang_rhys(delta_R: ndarray) float

Compute the Huang-Rhyes factor.

\[S_i = 1/2 \frac{\omega}{\hbar} [({M^{1/2}}^T \Delta R) \cdot \gamma_i]^2 = 1/2 \frac{\omega}{\hbar} {\sum_i m_i^{1/2} \gamma_i {\Delta R}_i}^2\]
Parameters:

delta_R – displacement in SI

localization_ratio()

Quantify the localization of the mode over the supercell.

See also participation_ratio().

1: fully delocalized >> 1: localized

participation_ratio()

Fraction of atoms active in the mode.

R J Bell et al 1970 J. Phys. C: Solid State Phys. 3 2111 https://doi.org/10.1088/0022-3719/3/10/013

It is equal to \(M_1^2 / (M_2 M_0)\) where

\[M_n = \sum_\alpha {m_\alpha {||\eta_\alpha||}^2}^n\]

where \(\eta_\alpha\) is the contribution of atom \(\alpha\) to eigendisplacement \(\eta\).

But \(M_0 = N\) by definition and \(M_1 = 1\) because the eigenvectors are normalized.

Note

\(M_n\) = np.sum(self.energies()**n)

per_species_n_eff(sp)

Compute the number of atoms participating to the mode for a given species.

per_species_pr(sp)

Compute the fraction of atoms participation to the mode for a given species.

See also per_species_n_eff().

project(delta_Q: ndarray) float

Project delta_Q onto the eigenvector.

project_coef2(delta_Q: ndarray) float

Square lenght of the projection of delta_Q onto the eigenvector.

project_coef2_R(delta_R: ndarray) float

Square lenght of the projection of delta_R onto the eigenvector.

set_lattice(lattice: ndarray, tol=1e-06) None

Change the representation to another lattice.

Parameters:
  • lattice – 3x3 matrix representing lattice vectors np.array([a, b, c]).

  • tol – numerical tolerance for vectors mismatch.

to_jmol(dest, **opts)

Write a mode into a Jmol file.

See hylight.jmol.export() for the parameters.

to_traj(duration, amplitude, framerate=25)

Produce a ase trajectory for animation purpose.

Parameters:
  • duration – duration of the animation in seconds

  • amplitude – amplitude applied to the mode in A (the modes are normalized)

  • framerate – number of frame per second of animation

hylight.mode.angle(v1, v2)

Compute the angle in radians between two 3D vectors.

hylight.mode.dynamical_matrix(phonons: Iterable[Mode]) ndarray

Retrieve the dynamical matrix from a set of modes.

Note that if some modes are missing the computation will fail.

Parameters:

phonons – list of modes

hylight.mode.generate_basis(seed)

Generate an orthonormal basis with the rows of seed as first rows.

Parameters:

seed – the starting vectors, a (m, n) matrix of orthonormal rows. m = 0 is valid and will create a random basis.

Returns:

a (n, n) orthonormal basis where the first m rows are the rows of seed.

hylight.mode.get_HR_factors(phonons: Iterable[Mode], delta_R_tot: ndarray, mask: Mask | None = None) ndarray

Compute the Huang-Rhys factors for all the real modes with energy above bias.

Parameters:
  • phonons – list of modes

  • delta_R_tot – displacement in SI

  • mask – a mask to filter modes based on their energies.

hylight.mode.get_energies(phonons: Iterable[Mode], mask: Mask | None = None) ndarray

Return an array of mode energies in SI.

hylight.mode.modes_from_dynmat(lattice, atoms, masses, ref, dynmat)

Compute vibrational modes from the dynamical matrix.

Parameters:
  • lattice – lattice parameters (3 x 3 numpy.ndarray)

  • atoms – list of atoms

  • masses – list of atom masses

  • ref – reference position

  • dynmat – dynamical matrix

Returns:

list of hylight.mode.Mode

hylight.mode.orthonormalize(m, n_skip=0)

Ensure that the vectors of m are orthonormal.

Change the rows from n_seed up inplace to make them orthonormal.

Parameters:
  • m – the starting vectors

  • n_seed – number of first rows to not change. They must be orthonormal already.

hylight.mode.project_on_asr(mat, masses)

Enforce accoustic sum rule.

Project the input matrix on the space of ASR abiding matrices and return the projections.

hylight.mode.rot_c_to_v(phonons: Iterable[Mode]) ndarray

Rotation matrix from Cartesian basis to Vibrational basis (right side).

hylight.mode.same_cell(cell1: ndarray, cell2: ndarray, tol=1e-06) CellMismatch | bool

Compare two lattice vectors matrix and return True if they describe the same cell.

hylight.mono_mode module

Simulation of luminescence spectra in 1D model.

hylight.mono_mode.compute_spectrum(e_zpl, S, sig, e_phonon_g, e=None)

Compute a spectrum from 1D model with experimental like inputs.

Parameters:
  • e_zpl – energy of the zero phonon line, in eV

  • S – the emission Huang-Rhys factor

  • sig – the lineshape standard deviation

  • e – (optional, None) a numpy array of energies to compute the spectrum at if ommited, an array ranging from 0 to 3*e_zpl will be created.

E_phonon_g:

the ground state phonon energy

hylight.mono_mode.huang_rhys(stokes_shift, e_phonon)

Huang-Rhys factor from the Stokes shift and the phonon energy.

hylight.multi_modes module

Simulation of spectra in nD model.

class hylight.multi_modes.LineShape(value, names=<not given>, *values, module=None, qualname=None, type=None, start=1, boundary=None)

Bases: Enum

Line shape type.

GAUSSIAN = 0
LORENTZIAN = 1
NONE = 2
hylight.multi_modes.compute_delta_R(poscar_gs, poscar_es)

Return \(\Delta R\) in A.

Parameters:
  • poscar_gs – path to ground state positions file.

  • poscar_es – path to excited state positions file.

Returns:

a numpy.ndarray of shape (n, 3) where n is the number of atoms.

hylight.multi_modes.compute_spectrum(phonons, delta_R, zpl, fwhm, e_max=None, resolution_e=0.001, mask=None, shape=LineShape.GAUSSIAN, pre_convolve=None, load_phonons=<function load_phonons>, window_fn=<function hamming>, T=0)

Compute a luminescence spectrum with the time-dependant formulation with an arbitrary linewidth.

Parameters:
  • phonons – list of modes (see load_phonons()) or path to load the modes

  • delta_R – displacement in A in a numpy array (see compute_delta_R()) or tuple of two paths (pos_gs, pos_es)

  • zpl – zero phonon line energy in eV

  • fwhm

    ZPL lineshape full width at half maximum in eV or None

    • if fwhm is None or fwhm == 0.0: the raw spectrum is provided unconvoluted

    • if fwhm > 0: the spectrum is convoluted with a gaussian line shape

    • if fwhm < 0: error

  • e_max – (optional, 2.5*e_zpl) max energy in eV (should be at greater than 2*zpl)

  • resolution_e – (optional, 1e-3) energy resolution in eV

  • load_phonons – a function to read phonons from files.

  • mask – a mode.Mask instance to select modes base on frequencies

  • shape – the lineshape (a LineShape instance)

  • pre_convolve – (float, optional, None) if not None, standard deviation of a pre convolution gaussian

  • window_fn – windowing function in the form provided by numpy (see numpy.hamming())

Returns:

(energy_array, intensity_array)

hylight.multi_modes.duschinsky(phonons_a, phonons_b)

Dushinsky matrix from b to a \(S_{a \\gets b}\).

hylight.multi_modes.dynmatshow(dynmat, blocks=None)

Plot the dynamical matrix.

Parameters:
  • dynmat – numpy array representing the dynamical matrice in SI.

  • blocks – (optional, None) a list of coloured blocks in the form (label, number_of_atoms, color).

hylight.multi_modes.fc_spectrum(phonons, delta_R, n_points=5000, disp=1)

Build arrays for plotting a spectrum energy spectral function.

hylight.multi_modes.freq_from_finite_diff(left, mid, right, mu, A=0.01)

Compute a vibration energy from three energy points.

Parameters:
  • left – energy (eV) of the left point

  • mid – energy (eV) of the middle point

  • right – energy (eV) of the right point

  • mu – effective mass associated with the displacement from the middle point to the sides.

  • A – amplitude (A) of the displacement between the middle point and the sides.

hylight.multi_modes.hr_spectrum(phonons, delta_R, n_points=5000, disp=1)

Build arrays for plotting a spectrum phonon spectral function.

hylight.multi_modes.make_line_shape(t, sigma_si, shape)

Create the lineshape function in time space.

Parameters:
  • t – the time array (in s)

  • sigma_si – the standard deviation of the line

  • shape – the type of lineshape (an instance of LineShape)

Returns:

a numpy.ndarray of the same shape as t

hylight.multi_modes.plot_spectral_function(mode_source, poscar_gs, poscar_es, load_phonons=<function load_phonons>, use_cm1=False, disp=1, mpl_params=None, mask=None)

Plot a two panel representation of the spectral function of the distorsion.

Parameters:
  • mode_source – path to the mode file (by default a pickle file)

  • poscar_gs – path to the file containing the ground state atomic positions.

  • poscar_es – path to the file containing the excited state atomic positions.

  • load_phonons – a function to read mode_source.

  • use_cm1 – use cm1 as the unit for frequency instead of meV.

  • disp – standard deviation of the gaussians in background in meV.

  • mpl_params – dictionary of kw parameters for pyplot.plot.

  • mask – a mask used in other computations to show on the plot.

Returns:

(figure, (ax_FC, ax_S))

hylight.multi_modes.rect(n)

A dummy windowing function that works like numpy.hamming, but as no effect on data.

hylight.npz module

Serialization of modes to numpy zipped file.

See also numpy.savez() and numpy.load().

hylight.npz.archive_modes(modes, dest, compress=False)

Store modes in dest using numpy’s npz format.

Parameters:
  • modes – a list of Mode objects.

  • dest – the path to write the modes to.

Returns:

the data returned by load_phonons.

hylight.npz.load_phonons(source)

Load modes from a Hylight archive.

hylight.npz.pops_and_masses(modes)

hylight.struct module

A generic representation of a crystal cell.

class hylight.struct.Struct(lattice, species, species_names=None)

Bases: object

A general description of a periodic crystal cell.

Store all the required infos to describe a given set of atomic positions.

property atoms

List the species names in an order matching self.raw.

copy() Struct

Return a copy of the structure.

classmethod from_mode(mode: Mode) Struct

Extract the cell from a hylight.mode.Mode.

get_offset(sp)

Return the index offset of the given species block in raw.

property raw

Return an array of atomic positions.

This can be modified overwritten, but not modified in place.

sp_at(i)

Return the name of the species at the given index.

Parameters:

i – the index of the atom.

Returns:

the name of the species.

property species_names

Names of species in the same order as found in the raw positions.

property system_name

The name of the system, eventually generated from formula.

Can be overwritten.

hylight.typing module

Helper module for type annotations.

Should help with variable support of typing across versions of python and libraries.

hylight.utils module

Pervasive utilities.

exception hylight.utils.InputError

Bases: ValueError

An exception raised when the input files are not as expected.

hylight.utils.gaussian(e, sigma, standard=True)

Evaluate a Gaussian function on e.

Parameters:
  • e – abscissa

  • sigma – standard deviation

  • standard

    • if True the curve is normalized to have an area of 1

    • if False the curve is normalized to have a maximum of 1

hylight.utils.gen_translat(lattice: ndarray)

Generate all translations to adjacent cells

Parameters:

lattice – np.ndarray([a, b, c]) first lattice parameter

hylight.utils.measure_fwhm(x, y)

Measure the full width at half maximum of a given spectrum.

Warning

It may fail if there are more than one band that reach half maximum in the array. In this case you may want to use select_interval to make a window around a single band.

Parameters:
  • x – the energy array

  • y – the intensity array

Returns:

FWHM in the same unit as x.

hylight.utils.parse_formatted_table(lines, format)

Parse a table of numbers from a list of string with a regex.

Parameters:
  • lines – a list of string representing the lines of the table

  • format – a re.Pattern or a str representing the format of the line It should fullmatch each line or a ValueError exception will be raised. All the groups defined in format will be converted to float64 by numpy.

Returns:

a np.array of dimension (len(lines), {number_of_groups})

Example:

>>> parse_formatted_table(["a=0.56 b=0.8 c=0.9"], "a=(.*) b=(.*) c=(.*)")
np.array([[0.56, 0.8, 0.9]])
hylight.utils.periodic_diff(lattice, ref, disp)

Compute the displacement between ref and disp, accounting for periodic conditions.

hylight.utils.periodic_dist(lattice, ref, disp)

Compute the distance between ref and disp, accounting for periodic conditions.

hylight.utils.select_interval(x, y, emin, emax, normalize=False, npoints=None)

Extract an interval of a spectrum and return the windows x and y arrays.

Parameters:
  • x – x array

  • y – y array

  • emin – lower bound for the window

  • emax – higher bound for the window

  • normalize – if true, the result y array is normalized

  • npoints – if an integer, the result arrays will be interpolated to contains exactly npoints linearly distributed between emin and emax.

Returns:

(windowed_x, windowed_y)

Module contents

The Hylight package.

Copyright (c) 2024, Théo Cavignac <theo.cavignac+dev@gmail.com>, The PyDEF team <camille.latouche@cnrs-imn.fr> Licensed under the EUPL

hylight.setup_logging()

Setup module logging.