pymatgen.phonon package
pymatgen phonon package with functionality for phonon DOS + bandstructure analysis and more.
Submodules
pymatgen.phonon.bandstructure module
This module provides classes to define a phonon band structure.
- class PhononBandStructure(qpoints: Sequence[Kpoint], frequencies: ArrayLike, lattice: Lattice, nac_frequencies: Sequence[Sequence] | None = None, eigendisplacements: ArrayLike = None, nac_eigendisplacements: Sequence[Sequence] | None = None, labels_dict: dict | None = None, coords_are_cartesian: bool = False, structure: Structure | None = None)[source]
Bases:
MSONable
This is the most generic phonon band structure data possible it’s defined by a list of qpoints + frequencies for each of them. Additional information may be given for frequencies at Gamma, where non-analytical contribution may be taken into account.
- Parameters:
qpoints – list of qpoint as numpy arrays, in frac_coords of the given lattice by default
frequencies – list of phonon frequencies in THz as a numpy array with shape (3*len(structure), len(qpoints)). The First index of the array refers to the band and the second to the index of the qpoint.
lattice – The reciprocal lattice as a pymatgen Lattice object. Pymatgen uses the physics convention of reciprocal lattice vectors WITH a 2*pi coefficient.
nac_frequencies – Frequencies with non-analytical contributions at Gamma in THz. A list of tuples. The first element of each tuple should be a list defining the direction (not necessarily a versor, will be normalized internally). The second element containing the 3*len(structure) phonon frequencies with non-analytical correction for that direction.
eigendisplacements – the phonon eigendisplacements associated to the frequencies in Cartesian coordinates. A numpy array of complex numbers with shape (3*len(structure), len(qpoints), len(structure), 3). The first index of the array refers to the band, the second to the index of the qpoint, the third to the atom in the structure and the fourth to the Cartesian coordinates.
nac_eigendisplacements – the phonon eigendisplacements associated to the non-analytical frequencies in nac_frequencies in Cartesian coordinates. A list of tuples. The first element of each tuple should be a list defining the direction. The second element containing a numpy array of complex numbers with shape (3*len(structure), len(structure), 3).
labels_dict – (dict[str, Kpoint]): this links a qpoint (in frac coords or Cartesian coordinates depending on the coords) to a label.
coords_are_cartesian (bool) – Whether the qpoint coordinates are Cartesian. Defaults to False.
structure – The crystal structure (as a pymatgen Structure object) associated with the band structure. This is needed to calculate element/orbital projections of the band structure.
- asr_breaking(tol_eigendisplacements: float = 1e-05) ndarray | None [source]
Get the breaking of the acoustic sum rule for the three acoustic modes, if Gamma is present. None otherwise. If eigendisplacements are available they are used to determine the acoustic modes: selects the bands corresponding to the eigendisplacements that represent to a translation within tol_eigendisplacements. If these are not identified or eigendisplacements are missing the first 3 modes will be used (indices [:3]).
- classmethod from_dict(dct: dict[str, Any]) Self [source]
- Parameters:
dct (dict) – Dict representation of PhononBandStructure.
- Returns:
PhononBandStructure
- get_gamma_point() Kpoint | None [source]
Get the Gamma q-point as a Kpoint object (or None if not found).
- get_nac_eigendisplacements_along_dir(direction) ndarray | None [source]
Get the nac_eigendisplacements for the given direction (not necessarily a versor). None if the direction is not present or nac_eigendisplacements has not been calculated.
- Parameters:
direction – the direction as a list of 3 elements
- Returns:
the eigendisplacements as a numpy array of complex numbers with shape (3*len(structure), len(structure), 3). None if not found.
- get_nac_frequencies_along_dir(direction: Sequence) np.ndarray | None [source]
Get the nac_frequencies for the given direction (not necessarily a versor). None if the direction is not present or nac_frequencies has not been calculated.
- Parameters:
direction – the direction as a list of 3 elements
- Returns:
the frequencies as a numpy array o(3*len(structure), len(qpoints)). None if not found.
- has_imaginary_freq(tol: float = 0.01) bool [source]
True if imaginary frequencies are present anywhere in the band structure. Always True if has_imaginary_gamma_freq is True.
- Parameters:
tol – Tolerance for determining if a frequency is imaginary. Defaults to 0.01.
- has_imaginary_gamma_freq(tol: float = 0.01) bool [source]
Check if there are imaginary modes at the gamma point and all close points.
- Parameters:
tol – Tolerance for determining if a frequency is imaginary. Defaults to 0.01.
- property has_nac: bool[source]
True if nac_frequencies are present (i.e. the band structure has been calculated taking into account Born-charge-derived non-analytical corrections at Gamma).
- max_freq() tuple[Kpoint, float] [source]
Get the q-point where the maximum frequency is reached and its value.
- class PhononBandStructureSymmLine(qpoints: Sequence[Kpoint], frequencies: ArrayLike, lattice: Lattice, has_nac: bool = False, eigendisplacements: ArrayLike = None, labels_dict: dict | None = None, coords_are_cartesian: bool = False, structure: Structure | None = None)[source]
Bases:
PhononBandStructure
Store phonon band structures along selected (symmetry) lines in the Brillouin zone. We call the different symmetry lines (ex: \Gamma to Z) “branches”.
- Parameters:
qpoints – list of qpoints as numpy arrays, in frac_coords of the given lattice by default
frequencies – list of phonon frequencies in eV as a numpy array with shape (3*len(structure), len(qpoints))
lattice – The reciprocal lattice as a pymatgen Lattice object. Pymatgen uses the physics convention of reciprocal lattice vectors WITH a 2*pi coefficient
has_nac – specify if the band structure has been produced taking into account non-analytical corrections at Gamma. If True frequencies at Gamma from different directions will be stored in naf. Default False.
eigendisplacements – the phonon eigendisplacements associated to the frequencies in Cartesian coordinates. A numpy array of complex numbers with shape (3*len(structure), len(qpoints), len(structure), 3). he First index of the array refers to the band, the second to the index of the qpoint, the third to the atom in the structure and the fourth to the Cartesian coordinates.
labels_dict – (dict) of {} this links a qpoint (in frac coords or Cartesian coordinates depending on the coords) to a label.
coords_are_cartesian – Whether the qpoint coordinates are cartesian.
structure – The crystal structure (as a pymatgen Structure object) associated with the band structure. This is needed if we provide projections to the band structure.
- as_phononwebsite() dict [source]
Return a dictionary with the phononwebsite format: https://henriquemiranda.github.io/phononwebsite.
- band_reorder() None [source]
Re-order the eigenvalues according to the similarity of the eigenvectors.
- classmethod from_dict(dct: dict) Self [source]
- Parameters:
dct (dict) – Dict representation.
- Returns:
PhononBandStructureSymmLine
- get_branch(index: int) list[dict[str, str | int]] [source]
Get in what branch(es) is the qpoint. There can be several branches.
- Parameters:
index (int) – the qpoint index
- Returns:
- [{“name”,”start_index”,”end_index”,”index”}]
indicating all branches in which the qpoint is. It takes into account the fact that one qpoint (e.g., \Gamma) can be in several branches
- Return type:
list[dict[str, str | int]]
- get_equivalent_qpoints(index: int) list[int] [source]
Get the list of qpoint indices equivalent (meaning they are the same frac coords) to the given one.
- Parameters:
index (int) – the qpoint index
- Returns:
equivalent indices
- Return type:
list[int]
TODO: now it uses the label we might want to use coordinates instead (in case there was a mislabel)
- write_phononwebsite(filename: str | PathLike) None [source]
Write a JSON file for the phononwebsite: https://henriquemiranda.github.io/phononwebsite.
- eigenvectors_from_displacements(disp: ndarray, masses: ndarray) ndarray [source]
Calculate the eigenvectors from the atomic displacements.
pymatgen.phonon.dos module
This module defines classes to represent the phonon density of states.
- class CompletePhononDos(structure: Structure, total_dos, ph_doses: dict)[source]
Bases:
PhononDos
This wrapper class defines a total dos, and also provides a list of PDos.
- pdos[source]
Dict of partial densities of the form {Site:Densities}. Densities are a dict of {Orbital:Values} where Values are a list of floats. Site is a pymatgen.core.sites.Site object.
- Type:
dict
- Parameters:
structure – Structure associated with this particular DOS.
total_dos – total Dos for structure
ph_doses – The phonon DOSes are supplied as a dict of {Site: Densities}.
- class PhononDos(frequencies: Sequence, densities: Sequence)[source]
Bases:
MSONable
Basic DOS object. All other DOS objects are extended versions of this object.
- Parameters:
frequencies – A sequence of frequencies in THz
densities – A sequence representing the density of states.
- cv(temp: float | None = None, structure: Structure | None = None, **kwargs) float [source]
Constant volume specific heat C_v at temperature T obtained from the integration of the DOS. Only positive frequencies will be used. Result in J/(K*mol-c). A mol-c is the abbreviation of a mole-cell, that is, the number of Avogadro times the atoms in a unit cell. To compare with experimental data the result should be divided by the number of unit formulas in the cell. If the structure is provided the division is performed internally and the result is in J/(K*mol).
- Parameters:
temp – a temperature in K
structure – the structure of the system. If not None it will be used to determine the number of formula units
**kwargs – allows passing in deprecated t parameter for temp
- Returns:
Constant volume specific heat C_v
- Return type:
float
- entropy(temp: float | None = None, structure: Structure | None = None, **kwargs) float [source]
Vibrational entropy at temperature T obtained from the integration of the DOS. Only positive frequencies will be used. Result in J/(K*mol-c). A mol-c is the abbreviation of a mole-cell, that is, the number of Avogadro times the atoms in a unit cell. To compare with experimental data the result should be divided by the number of unit formulas in the cell. If the structure is provided the division is performed internally and the result is in J/(K*mol).
- Parameters:
temp – a temperature in K
structure – the structure of the system. If not None it will be used to determine the number of formula units
**kwargs – allows passing in deprecated t parameter for temp
- Returns:
Vibrational entropy
- Return type:
float
- static fp_to_dict(fp: PhononDosFingerprint) dict [source]
Convert a fingerprint into a dictionary.
- Parameters:
fp – The DOS fingerprint to be converted into a dictionary
- Returns:
A dict of the fingerprint Keys=type, Values=np.ndarray(frequencies, densities)
- Return type:
dict
- classmethod from_dict(dct: dict[str, Sequence]) Self [source]
Get PhononDos object from dict representation of PhononDos.
- get_dos_fp(binning: bool = True, min_f: float | None = None, max_f: float | None = None, n_bins: int = 256, normalize: bool = True) PhononDosFingerprint [source]
Generate the DOS fingerprint.
- Parameters:
binning (bool) – If true, the DOS fingerprint is binned using np.linspace and n_bins. Default is True.
min_f (float) – The minimum mode frequency to include in the fingerprint (default is None)
max_f (float) – The maximum mode frequency to include in the fingerprint (default is None)
n_bins (int) – Number of bins to be used in the fingerprint (default is 256)
normalize (bool) – If true, normalizes the area under fp to equal to 1. Default is True.
- Returns:
- The phonon density of states fingerprint
of format (frequencies, densities, n_bins, bin_width)
- Return type:
- static get_dos_fp_similarity(fp1: PhononDosFingerprint, fp2: PhononDosFingerprint, col: int = 1, pt: int | str = 'All', normalize: bool = False, metric: Literal['tanimoto', 'wasserstein', 'cosine-sim'] = 'tanimoto') float [source]
Calculate the similarity index between two fingerprints.
- Parameters:
fp1 (PhononDosFingerprint) – The 1st dos fingerprint object
fp2 (PhononDosFingerprint) – The 2nd dos fingerprint object
col (int) – The item in the fingerprints (0:frequencies,1: densities) to compute the similarity index of (default is 1)
pt (int or str) – The index of the point that the dot product is to be taken (default is All)
normalize (bool) – If True normalize the scalar product to 1 (default is False)
metric (Literal) – Metric used to compute similarity default is “tanimoto”.
- Raises:
ValueError – If metric other than “tanimoto”, “wasserstein” and “cosine-sim” is requested.
ValueError – If normalize is set to True along with the metric.
- Returns:
Similarity index given by the dot product
- Return type:
float
- get_interpolated_value(frequency) float [source]
Get interpolated density for a particular frequency.
- Parameters:
frequency – frequency to return the density for.
- get_last_peak(threshold: float = 0.05) float [source]
Find the last peak in the phonon DOS defined as the highest frequency with a DOS value at least threshold * height of the overall highest DOS peak. A peak is any local maximum of the DOS as a function of frequency. Use dos.get_interpolated_value(peak_freq) to get density at peak_freq.
TODO method added by @janosh on 2023-12-18. seems to work in most cases but was not extensively tested. PRs with improvements welcome!
- Parameters:
threshold (float, optional) – Minimum ratio of the height of the last peak to the height of the highest peak. Defaults to 0.05 = 5%. In case no peaks are high enough to match, the threshold is reset to half the height of the second-highest peak.
- Returns:
last DOS peak frequency (in THz)
- Return type:
float
- get_smeared_densities(sigma: float) ndarray [source]
Get the densities, but with a Gaussian smearing of std dev sigma applied.
- Parameters:
sigma – Std dev of Gaussian smearing function. In units of THz. Common values are 0.01 - 0.1 THz.
- Returns:
Gaussian-smeared DOS densities.
- Return type:
np.array
- helmholtz_free_energy(temp: float | None = None, structure: Structure | None = None, **kwargs) float [source]
Phonon contribution to the Helmholtz free energy at temperature T obtained from the integration of the DOS. Only positive frequencies will be used. Result in J/mol-c. A mol-c is the abbreviation of a mole-cell, that is, the number of Avogadro times the atoms in a unit cell. To compare with experimental data the result should be divided by the number of unit formulas in the cell. If the structure is provided the division is performed internally and the result is in J/mol.
- Parameters:
temp – a temperature in K
structure – the structure of the system. If not None it will be used to determine the number of formula units
**kwargs – allows passing in deprecated t parameter for temp
- Returns:
Phonon contribution to the Helmholtz free energy
- Return type:
float
- internal_energy(temp: float | None = None, structure: Structure | None = None, **kwargs) float [source]
Phonon contribution to the internal energy at temperature T obtained from the integration of the DOS. Only positive frequencies will be used. Result in J/mol-c. A mol-c is the abbreviation of a mole-cell, that is, the number of Avogadro times the atoms in a unit cell. To compare with experimental data the result should be divided by the number of unit formulas in the cell. If the structure is provided the division is performed internally and the result is in J/mol.
- Parameters:
temp – a temperature in K
structure – the structure of the system. If not None it will be used to determine the number of formula units
**kwargs – allows passing in deprecated t parameter for temp
- Returns:
Phonon contribution to the internal energy
- Return type:
float
- mae(other: PhononDos, two_sided: bool = True) float [source]
Mean absolute error between two DOSs.
- Parameters:
other (PhononDos) – Another phonon DOS
two_sided (bool) – Whether to calculate the two-sided MAE meaning interpolate each DOS to the other’s frequencies and averaging the two MAEs. Defaults to True.
- Returns:
Mean absolute error.
- Return type:
float
- r2_score(other: PhononDos) float [source]
R^2 score between two DOSs.
- Parameters:
other (PhononDos) – Another phonon DOS
- Returns:
R^2 score
- Return type:
float
- zero_point_energy(structure: Structure | None = None) float [source]
Zero point energy of the system. Only positive frequencies will be used. Result in J/mol-c. A mol-c is the abbreviation of a mole-cell, that is, the number of Avogadro times the atoms in a unit cell. To compare with experimental data the result should be divided by the number of unit formulas in the cell. If the structure is provided the division is performed internally and the result is in J/mol.
- Parameters:
structure – the structure of the system. If not None it will be used to determine the number of formula units
- Returns:
Phonon contribution to the internal energy
- class PhononDosFingerprint(frequencies: NDArray, densities: NDArray, n_bins: int, bin_width: float)[source]
Bases:
NamedTuple
Represents a Phonon Density of States (DOS) fingerprint.
This named tuple is used to store information related to the Density of States (DOS) in a material. It includes the frequencies, densities, number of bins, and bin width.
- Parameters:
frequencies – The frequency values associated with the DOS.
densities – The corresponding density values for each energy.
n_bins – The number of bins used in the fingerprint.
bin_width – The width of each bin in the DOS fingerprint.
Create new instance of PhononDosFingerprint(frequencies, densities, n_bins, bin_width)
pymatgen.phonon.gruneisen module
This module provides classes to define a Grueneisen band structure.
- class GruneisenParameter(qpoints: ArrayLike, gruneisen: ArrayLike[ArrayLike], frequencies: ArrayLike[ArrayLike], multiplicities: Sequence | None = None, structure: Structure = None, lattice: Lattice = None)[source]
Bases:
MSONable
Store the Gruneisen parameter for a single q-point on a regular grid.
- Parameters:
qpoints – list of qpoints as numpy arrays, in frac_coords of the given lattice by default
gruneisen – list of gruneisen parameters as numpy arrays, shape: (3*len(structure), len(qpoints))
frequencies – list of phonon frequencies in THz as a numpy array with shape (3*len(structure), len(qpoints))
multiplicities – list of multiplicities
structure – The crystal structure (as a pymatgen Structure object) associated with the gruneisen parameters.
lattice – The reciprocal lattice as a pymatgen Lattice object. Pymatgen uses the physics convention of reciprocal lattice vectors WITH a 2*pi coefficient.
- property acoustic_debye_temp: float[source]
Acoustic Debye temperature in K, i.e. the Debye temperature divided by n_sites**(1/3). Adapted from abipy.
- average_gruneisen(t: float | None = None, squared: bool = True, limit_frequencies: Literal['debye', 'acoustic'] | None = None) float [source]
Calculate the average of the Gruneisen based on the values on the regular grid. If squared is True, the average will use the squared value of the Gruneisen and a squared root is performed on the final result. Values associated with negative frequencies will be ignored. See Nath et al. _Scripta Materialia_ 2017, _129_, 88 for the definitions. Adapted from classes in abipy that have been written by Guido Petretto (UCLouvain).
- Parameters:
t – the temperature at which the average Gruneisen will be evaluated. If None the acoustic Debye temperature is used (see acoustic_debye_temp).
squared – if True the average is performed on the squared values of the Grueneisen.
limit_frequencies – if None (default) no limit on the frequencies will be applied. Possible values are “debye” (only modes with frequencies lower than the acoustic Debye temperature) and “acoustic” (only the acoustic modes, i.e. the first three modes).
- Returns:
The average Gruneisen parameter
- debye_temp_phonopy(freq_max_fit=None) float [source]
Get Debye temperature in K as implemented in phonopy.
- Parameters:
freq_max_fit – Maximum frequency to include for fitting. Defaults to include first quartile of frequencies.
- Returns:
Debye temperature in K.
- thermal_conductivity_slack(squared: bool = True, limit_frequencies: Literal['debye', 'acoustic'] | None = None, theta_d: float | None = None, t: float | None = None) float [source]
Calculate the thermal conductivity at the acoustic Debye temperature with the Slack formula, using the average Gruneisen. Adapted from abipy.
- Parameters:
squared (bool) – if True the average is performed on the squared values of the Gruneisen
limit_frequencies – if None (default) no limit on the frequencies will be applied. Possible values are “debye” (only modes with frequencies lower than the acoustic Debye temperature) and “acoustic” (only the acoustic modes, i.e. the first three modes).
theta_d – the temperature used to estimate the average of the Gruneisen used in the Slack formula. If None the acoustic Debye temperature is used (see acoustic_debye_temp). Will also be considered as the Debye temperature in the Slack formula.
t – temperature at which the thermal conductivity is estimated. If None the value at the calculated acoustic Debye temperature is given. The value is obtained as a simple rescaling of the value at the Debye temperature.
- Returns:
The value of the thermal conductivity in W/(m*K)
- class GruneisenPhononBandStructure(qpoints: ArrayLike, frequencies: ArrayLike[ArrayLike], gruneisenparameters: ArrayLike, lattice: Lattice, eigendisplacements: ArrayLike[ArrayLike] = None, labels_dict: dict | None = None, coords_are_cartesian: bool = False, structure: Structure | None = None)[source]
Bases:
PhononBandStructure
This is the most generic phonon band structure data possible it’s defined by a list of qpoints + frequencies for each of them. Additional information may be given for frequencies at Gamma, where non-analytical contribution may be taken into account.
- Parameters:
qpoints – list of qpoint as numpy arrays, in frac_coords of the given lattice by default
frequencies – list of phonon frequencies in THz as a numpy array with shape (3*len(structure), len(qpoints)). The First index of the array refers to the band and the second to the index of the qpoint.
gruneisenparameters – list of Grueneisen parameters with the same structure frequencies.
lattice – The reciprocal lattice as a pymatgen Lattice object. Pymatgen uses the physics convention of reciprocal lattice vectors WITH a 2*pi coefficient.
eigendisplacements – the phonon eigendisplacements associated to the frequencies in Cartesian coordinates. A numpy array of complex numbers with shape (3*len(structure), len(qpoints), len(structure), 3). The first index of the array refers to the band, the second to the index of the qpoint, the third to the atom in the structure and the fourth to the Cartesian coordinates.
labels_dict – (dict) of {} this links a qpoint (in frac coords or Cartesian coordinates depending on the coords) to a label.
coords_are_cartesian – Whether the qpoint coordinates are Cartesian.
structure – The crystal structure (as a pymatgen Structure object) associated with the band structure. This is needed if we provide projections to the band structure.
- class GruneisenPhononBandStructureSymmLine(qpoints: ArrayLike, frequencies: ArrayLike[ArrayLike], gruneisenparameters: ArrayLike, lattice: Lattice, eigendisplacements: ArrayLike[ArrayLike] = None, labels_dict: dict | None = None, coords_are_cartesian: bool = False, structure: Structure | None = None)[source]
Bases:
GruneisenPhononBandStructure
,PhononBandStructureSymmLine
Store a GruneisenPhononBandStructureSymmLine together with Grueneisen parameters for every frequency.
- Parameters:
qpoints – list of qpoints as numpy arrays, in frac_coords of the given lattice by default
frequencies – list of phonon frequencies in eV as a numpy array with shape (3*len(structure), len(qpoints))
gruneisenparameters – list of Grueneisen parameters as a numpy array with the shape (3*len(structure), len(qpoints))
lattice – The reciprocal lattice as a pymatgen Lattice object. Pymatgen uses the physics convention of reciprocal lattice vectors WITH a 2*pi coefficient
eigendisplacements – the phonon eigendisplacements associated to the frequencies in Cartesian coordinates. A numpy array of complex numbers with shape (3*len(structure), len(qpoints), len(structure), 3). The first index of the array refers to the band, the second to the index of the qpoint, the third to the atom in the structure and the fourth to the Cartesian coordinates.
labels_dict – (dict) of {} this links a qpoint (in frac coords or Cartesian coordinates depending on the coords) to a label.
coords_are_cartesian – Whether the qpoint coordinates are cartesian.
structure – The crystal structure (as a pymatgen Structure object) associated with the band structure. This is needed if we provide projections to the band structure.
pymatgen.phonon.ir_spectra module
This module provides classes to handle the calculation of the IR spectra This implementation is adapted from Abipy https://github.com/abinit/abipy where it was originally done by Guido Petretto and Matteo Giantomassi.
- class IRDielectricTensor(oscillator_strength: ArrayLike, ph_freqs_gamma: ArrayLike, epsilon_infinity: ArrayLike, structure: Structure)[source]
Bases:
MSONable
Handle the Ionic Dielectric Tensor The implementation is adapted from Abipy See the definitions Eq.(53-54) in :cite:`Gonze1997` PRB55, 10355 (1997).
- Parameters:
oscillator_strength – IR oscillator strengths as defined in Eq. 54 in :cite:`Gonze1997` PRB55, 10355 (1997).
ph_freqs_gamma – Phonon frequencies at the Gamma point
epsilon_infinity – electronic susceptibility as defined in Eq. 29.
structure – A Structure object corresponding to the structure used for the calculation.
- get_ir_spectra(broad: list | float = 5e-05, emin: float = 0, emax: float | None = None, divs: int = 500) tuple [source]
The IR spectra is obtained for the different directions.
- Parameters:
broad – a list of broadenings or a single broadening for the phonon peaks
emin (float) – minimum energy in which to obtain the spectra. Defaults to 0.
emax (float) – maximum energy in which to obtain the spectra. Defaults to None.
divs – number of frequency samples between emin and emax
- Returns:
- divs array with the frequencies at which the
dielectric tensor is calculated
- dielectric_tensor: divsx3x3 numpy array with the dielectric tensor
for the range of frequencies
- Return type:
frequencies
- get_plotter(components: Sequence = ('xx',), reim: str = 'reim', broad: list | float = 5e-05, emin: float = 0, emax: float | None = None, divs: int = 500, **kwargs) SpectrumPlotter [source]
Return an instance of the Spectrum plotter containing the different requested components.
- Parameters:
components – A list with the components of the dielectric tensor to plot. Can be either two indexes or a string like ‘xx’ to plot the (0,0) component
reim – If ‘re’ (im) is present in the string plots the real (imaginary) part of the dielectric tensor
broad (float) – a list of broadenings or a single broadening for the phonon peaks. Defaults to 0.00005.
emin (float) – minimum energy in which to obtain the spectra. Defaults to 0.
emax (float) – maximum energy in which to obtain the spectra. Defaults to None.
divs – number of frequency samples between emin and emax
**kwargs – Passed to IRDielectricTensor.get_spectrum()
- get_spectrum(component: Sequence | str, reim: str, broad: list | float = 5e-05, emin: float = 0, emax: float | None = None, divs: int = 500, label=None) Spectrum [source]
component: either two indexes or a string like ‘xx’ to plot the (0,0) component reim: only “re” or “im” broad: a list of broadenings or a single broadening for the phonon peaks.
- plot(components: Sequence = ('xx',), reim: str = 'reim', show_phonon_frequencies: bool = True, xlim: float | None = None, ylim: float | None = None, **kwargs) Axes [source]
Helper function to generate the Spectrum plotter and directly plot the results.
- Parameters:
components – A list with the components of the dielectric tensor to plot. Can be either two indexes or a string like ‘xx’ to plot the (0,0) component
reim – If ‘re’ (im) is present in the string plots the real (imaginary) part of the dielectric tensor
show_phonon_frequencies – plot a dot where the phonon frequencies are to help identify IR inactive modes
xlim – x-limits of the plot. Defaults to None for automatic determination.
ylim – y-limits of the plot. Defaults to None for automatic determination.
kwargs – keyword arguments passed to the plotter
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.
pymatgen.phonon.plotter module
This module implements plotter for DOS and band structure.
- class FreqUnits(factor: float, label: str)[source]
Bases:
NamedTuple
Conversion factor from THz to the required units and the label.
Create new instance of FreqUnits(factor, label)
- class GruneisenPhononBSPlotter(bs: GruneisenPhononBandStructureSymmLine)[source]
Bases:
PhononBSPlotter
Plot or get data to facilitate the plot of band structure objects.
- Parameters:
bs – A GruneisenPhononBandStructureSymmLine object.
- bs_plot_data() dict[str, Any] [source]
Get the data nicely formatted for a plot.
- Returns:
ticks: A dict with the ‘distances’ at which there is a qpoint (the x axis) and the labels (None if no label) frequencies: A list (one element for each branch) of frequencies for each qpoint: [branch][qpoint][mode]. The data is stored by branch to facilitate the plotting gruneisen: GruneisenPhononBandStructureSymmLine lattice: The reciprocal lattice.
- Return type:
A dict of the following format
- get_plot_gs(ylim: float | None = None, plot_ph_bs_with_gruneisen: bool = False, **kwargs) Axes [source]
Get a matplotlib object for the Gruneisen bandstructure plot.
- Parameters:
ylim – Specify the y-axis (gruneisen) limits; by default None let the code choose.
plot_ph_bs_with_gruneisen (bool) – Plot phonon band-structure with bands coloured as per Gruneisen parameter values on a logarithmic scale
**kwargs – additional keywords passed to ax.plot().
- plot_compare_gs(other_plotter: GruneisenPhononBSPlotter) Axes [source]
Plot two band structure for comparison. One is in red the other in blue. The two band structures need to be defined on the same symmetry lines! and the distance between symmetry lines is the one of the band structure used to build the PhononBSPlotter.
- Parameters:
other_plotter (GruneisenPhononBSPlotter) – another phonon DOS plotter defined along the same symmetry lines.
- Raises:
ValueError – if the two plotters are incompatible (due to different data lengths)
- Returns:
with both band structures
- Return type:
plt.Axes
- save_plot_gs(filename: str | PathLike, img_format: str = 'eps', ylim: float | None = None, plot_ph_bs_with_gruneisen: bool = False, **kwargs) None [source]
Save matplotlib plot to a file.
- Parameters:
filename – Filename to write to.
img_format – Image format to use. Defaults to EPS.
ylim – Specifies the y-axis limits.
plot_ph_bs_with_gruneisen – Plot phonon band-structure with bands coloured as per Gruneisen parameter values on a logarithmic scale
**kwargs – kwargs passed to get_plot_gs
- show_gs(ylim: float | None = None, plot_ph_bs_with_gruneisen: bool = False, **kwargs) None [source]
Show the plot using matplotlib.
- Parameters:
ylim – Specifies the y-axis limits.
plot_ph_bs_with_gruneisen – Plot phonon band-structure with bands coloured as per Gruneisen parameter values on a logarithmic scale
**kwargs – kwargs passed to get_plot_gs
- class GruneisenPlotter(gruneisen: GruneisenParameter)[source]
Bases:
object
Plot GruneisenParameter.
Plot information from GruneisenParameter.
- Parameters:
gruneisen (GruneisenParameter) – containing the data to plot.
- get_plot(marker: str = 'o', markersize: float | None = None, units: Literal['thz', 'ev', 'mev', 'ha', 'cm-1', 'cm^-1'] = 'thz') Axes [source]
Will produce a plot.
- Parameters:
marker – marker for the depiction
markersize – size of the marker
units – unit for the plots, accepted units: thz, ev, mev, ha, cm-1, cm^-1.
- Returns:
matplotlib axes object
- Return type:
plt.Axes
- save_plot(filename: str | PathLike, img_format: str = 'pdf', units: Literal['thz', 'ev', 'mev', 'ha', 'cm-1', 'cm^-1'] = 'thz') None [source]
Will save the plot to a file.
- Parameters:
filename – name of the filename
img_format – format of the saved plot
units – accepted units: thz, ev, mev, ha, cm-1, cm^-1.
- class PhononBSPlotter(bs: PhononBandStructureSymmLine, label: str | None = None)[source]
Bases:
object
Plot or get data to facilitate the plot of band structure objects.
- Parameters:
bs – A PhononBandStructureSymmLine object.
label – A label for the plot. Defaults to None for no label. Esp. useful with the plot_compare method to distinguish the band structures.
- bs_plot_data() dict[str, Any] [source]
Get the data nicely formatted for a plot.
- Returns:
ticks: A dict with the ‘distances’ at which there is a qpoint (the x axis) and the labels (None if no label) frequencies: A list (one element for each branch) of frequencies for each qpoint: [branch][qpoint][mode]. The data is stored by branch to facilitate the plotting lattice: The reciprocal lattice.
- Return type:
A dict of the following format
- get_plot(ylim: float | None = None, units: Literal['thz', 'ev', 'mev', 'ha', 'cm-1', 'cm^-1'] = 'thz', **kwargs) Axes [source]
Get a matplotlib object for the bandstructure plot.
- Parameters:
ylim – Specify the y-axis (frequency) limits; by default None let the code choose.
units – units for the frequencies. Accepted values thz, ev, mev, ha, cm-1, cm^-1. Defaults to “thz”.
**kwargs – passed to ax.plot function.
- get_proj_plot(site_comb: str | list[list[int]] = 'element', ylim: tuple[None | float, None | float] | None = None, units: Literal['thz', 'ev', 'mev', 'ha', 'cm-1', 'cm^-1'] = 'thz', rgb_labels: tuple[None | str] | None = None) Axes [source]
Get a matplotlib object for the bandstructure plot projected along atomic sites.
- Parameters:
site_comb – a list of list, for example, [[0],[1],[2,3,4]]; the numbers in each sublist represents the indices of atoms; the atoms in a same sublist will be plotted in a same color; if not specified, unique elements are automatically grouped.
ylim – Specify the y-axis (frequency) limits; by default None let the code choose.
units – units for the frequencies. Accepted values thz, ev, mev, ha, cm-1, cm^-1. Defaults to “thz”.
rgb_labels – a list of rgb colors for the labels; if not specified, the colors will be automatically generated.
- get_ticks() dict[str, list] [source]
Get all ticks and labels for a band structure plot.
- Returns:
a list of distance at which ticks should be set and ‘label’: a list of label for each of those ticks.
- Return type:
A dict with ‘distance’
- plot_compare(other_plotter: PhononBSPlotter | dict[str, PhononBSPlotter], units: Literal['thz', 'ev', 'mev', 'ha', 'cm-1', 'cm^-1'] = 'thz', self_label: str = 'self', colors: Sequence[str] | None = None, legend_kwargs: dict | None = None, on_incompatible: Literal['raise', 'warn', 'ignore'] = 'raise', other_kwargs: dict | None = None, **kwargs) Axes [source]
Plot two band structure for comparison. self in blue, others in red, green, … The band structures need to be defined on the same symmetry lines! The distance between symmetry lines is determined by the band structure used to initialize PhononBSPlotter (self).
- Parameters:
other_plotter (PhononBSPlotter | dict[str, PhononBSPlotter]) – Other PhononBSPlotter object(s) defined along the same symmetry lines
units (str) – units for the frequencies. Accepted values thz, ev, mev, ha, cm-1, cm^-1. Defaults to ‘thz’.
self_label (str) – label for the self band structure. Defaults to to the label passed to PhononBSPlotter.init or, if None, ‘self’.
colors (list[str]) – list of colors for the other band structures. Defaults to None for automatic colors.
legend_kwargs – dict[str, Any]: kwargs passed to ax.legend().
on_incompatible ('raise' | 'warn' | 'ignore') – What to do if the band structures are not compatible. Defaults to ‘raise’.
other_kwargs – dict[str, Any]: kwargs passed to other_plotter ax.plot().
**kwargs – passed to ax.plot().
- Returns:
with two band structures.
- Return type:
plt.Axes
- save_plot(filename: str | PathLike, ylim: float | None = None, units: Literal['thz', 'ev', 'mev', 'ha', 'cm-1', 'cm^-1'] = 'thz') None [source]
Save matplotlib plot to a file.
- Parameters:
filename (str | Path) – Filename to write to.
ylim (float) – Specifies the y-axis limits.
units ("thz" | "ev" | "mev" | "ha" | "cm-1" | "cm^-1") – units for the frequencies.
- show(ylim: float | None = None, units: Literal['thz', 'ev', 'mev', 'ha', 'cm-1', 'cm^-1'] = 'thz') None [source]
Show the plot using matplotlib.
- Parameters:
ylim (float) – Specifies the y-axis limits.
units ("thz" | "ev" | "mev" | "ha" | "cm-1" | "cm^-1") – units for the frequencies.
- show_proj(site_comb: str | list[list[int]] = 'element', ylim: tuple[None | float, None | float] | None = None, units: Literal['thz', 'ev', 'mev', 'ha', 'cm-1', 'cm^-1'] = 'thz', rgb_labels: tuple[str] | None = None) None [source]
Show the projected plot using matplotlib.
- Parameters:
site_comb – A list of list of indices of sites to combine. For example, [[0, 1], [2, 3]] will combine the projections of sites 0 and 1, and sites 2 and 3. Defaults to “element”, which will combine sites by element.
ylim – Specify the y-axis (frequency) limits; by default None let the code choose.
units – units for the frequencies. Accepted values thz, ev, mev, ha, cm-1, cm^-1. Defaults to “thz”.
rgb_labels – A list of labels for the rgb triangle. Defaults to None, which will use the element symbols.
- class PhononDosPlotter(stack: bool = False, sigma: float | None = None)[source]
Bases:
object
Plot phonon DOSs. The interface is very flexible to accommodate the many different ways in which people want to view DOS. Typical usage is:
# Initializes plotter with some optional args. Defaults are usually fine plotter = PhononDosPlotter().
# Add DOS with a label plotter.add_dos(“Total DOS”, dos)
# Alternatively, you can add a dict of DOSes. This is the typical form # returned by CompletePhononDos.get_element_dos().
- Parameters:
stack – Whether to plot the DOS as a stacked area graph
sigma – A float specifying a standard deviation for Gaussian smearing the DOS for nicer looking plots. Defaults to None for no smearing.
- add_dos(label: str, dos: PhononDos, **kwargs: Any) None [source]
Add a dos for plotting.
- Parameters:
label (str) – label for the DOS. Must be unique.
dos (PhononDos) – DOS object
**kwargs – kwargs supported by matplotlib.pyplot.plot
- add_dos_dict(dos_dict: dict, key_sort_func=None) None [source]
Add a dictionary of doses, with an optional sorting function for the keys.
- Parameters:
dos_dict – dict of {label: Dos}
key_sort_func – function used to sort the dos_dict keys.
- get_dos_dict() dict [source]
Get the added doses as a json-serializable dict. Note that if you have specified smearing for the DOS plot, the densities returned will be the smeared densities, not the original densities.
- Returns:
DOS data. Generally of the form {label: {‘frequencies’:.., ‘densities’: …}}
- Return type:
dict
- get_plot(xlim: float | None = None, ylim: float | None = None, invert_axes: bool = False, units: Literal['thz', 'ev', 'mev', 'ha', 'cm-1', 'cm^-1'] = 'thz', legend: dict | None = None, ax: Axes | None = None) Axes [source]
Get a matplotlib plot showing the DOS.
- Parameters:
xlim – Specifies the x-axis limits. Set to None for automatic determination.
ylim – Specifies the y-axis limits.
invert_axes (bool) – Whether to invert the x and y axes. Enables chemist style DOS plotting. Defaults to False.
units (thz | ev | mev | ha | cm-1 | cm^-1) – units for the frequencies. Defaults to “thz”.
legend – dict with legend options. For example, {“loc”: “upper right”} will place the legend in the upper right corner. Defaults to {“fontsize”: 30}.
ax (Axes) – An existing axes object onto which the plot will be added. If None, a new figure will be created.
- save_plot(filename: str | PathLike, img_format: str = 'eps', xlim: float | None = None, ylim: float | None = None, invert_axes: bool = False, units: Literal['thz', 'ev', 'mev', 'ha', 'cm-1', 'cm^-1'] = 'thz') None [source]
Save matplotlib plot to a file.
- Parameters:
filename – Filename to write to.
img_format – Image format to use. Defaults to EPS.
xlim – Specifies the x-axis limits. Set to None for automatic determination.
ylim – Specifies the y-axis limits.
invert_axes – Whether to invert the x and y axes. Enables chemist style DOS plotting. Defaults to False.
units – units for the frequencies. Accepted values thz, ev, mev, ha, cm-1, cm^-1
- show(xlim: float | None = None, ylim: None = None, invert_axes: bool = False, units: Literal['thz', 'ev', 'mev', 'ha', 'cm-1', 'cm^-1'] = 'thz') None [source]
Show the plot using matplotlib.
- Parameters:
xlim – Specifies the x-axis limits. Set to None for automatic determination.
ylim – Specifies the y-axis limits.
invert_axes – Whether to invert the x and y axes. Enables chemist style DOS plotting. Defaults to False.
units – units for the frequencies. Accepted values thz, ev, mev, ha, cm-1, cm^-1.
- class ThermoPlotter(dos: PhononDos, structure: Structure = None)[source]
Bases:
object
Plotter for thermodynamic properties obtained from phonon DOS. If the structure corresponding to the DOS, it will be used to extract the formula unit and provide the plots in units of mol instead of mole-cell.
- Parameters:
dos – A PhononDos object.
structure – A Structure object corresponding to the structure used for the calculation.
- plot_cv(tmin: float, tmax: float, ntemp: int, ylim: float | None = None, **kwargs) Figure [source]
Plots the constant volume specific heat C_v in a temperature range.
- Parameters:
tmin – minimum temperature
tmax – maximum temperature
ntemp – number of steps
ylim – tuple specifying the y-axis limits.
kwargs – kwargs passed to the matplotlib function ‘plot’.
- Returns:
matplotlib figure
- Return type:
plt.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.
- plot_entropy(tmin: float, tmax: float, ntemp: int, ylim: float | None = None, **kwargs) Figure [source]
Plots the vibrational entrpy in a temperature range.
- Parameters:
tmin – minimum temperature
tmax – maximum temperature
ntemp – number of steps
ylim – tuple specifying the y-axis limits.
kwargs – kwargs passed to the matplotlib function ‘plot’.
- Returns:
matplotlib figure
- Return type:
plt.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.
- plot_helmholtz_free_energy(tmin: float, tmax: float, ntemp: int, ylim: float | None = None, **kwargs) Figure [source]
Plots the vibrational contribution to the Helmoltz free energy in a temperature range.
- Parameters:
tmin – minimum temperature
tmax – maximum temperature
ntemp – number of steps
ylim – tuple specifying the y-axis limits.
kwargs – kwargs passed to the matplotlib function ‘plot’.
- Returns:
matplotlib figure
- Return type:
plt.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.
- plot_internal_energy(tmin: float, tmax: float, ntemp: int, ylim: float | None = None, **kwargs) Figure [source]
Plots the vibrational internal energy in a temperature range.
- Parameters:
tmin – minimum temperature
tmax – maximum temperature
ntemp – number of steps
ylim – tuple specifying the y-axis limits.
kwargs – kwargs passed to the matplotlib function ‘plot’.
- Returns:
matplotlib figure
- Return type:
plt.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.
- plot_thermodynamic_properties(tmin: float, tmax: float, ntemp: int, ylim: float | None = None, **kwargs) Figure [source]
Plots all the thermodynamic properties in a temperature range.
- Parameters:
tmin – minimum temperature
tmax – maximum temperature
ntemp – number of steps
ylim – tuple specifying the y-axis limits.
kwargs – kwargs passed to the matplotlib function ‘plot’.
- Returns:
matplotlib figure
- Return type:
plt.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.
pymatgen.phonon.thermal_displacements module
This module provides classes to handle thermal displacement matrices (anisotropic displacement parameters).
- class ThermalDisplacementMatrices(thermal_displacement_matrix_cart: ArrayLike[ArrayLike], structure: Structure, temperature: float | None, thermal_displacement_matrix_cif: ArrayLike[ArrayLike] = None)[source]
Bases:
MSONable
Handle thermal displacement matrices This class stores thermal displacement matrices in Ucart format.
An earlier implementation based on Matlab can be found here: https://github.com/JaGeo/MolecularToolbox ( J. George, A. Wang, V. L. Deringer, R. Wang, R. Dronskowski, U. Englert, CrystEngComm, 2015, 17, 7414-7422.)
- Parameters:
thermal_displacement_matrix_cart – 2D numpy array including the thermal_displacement matrix Ucart 1st dimension atom types, then compressed thermal displacement matrix will follow U11, U22, U33, U23, U13, U12 (xx, yy, zz, yz, xz, xy) convention similar to “thermal_displacement_matrices.yaml” in phonopy
structure – A pymatgen Structure object
temperature – temperature at which thermal displacement matrix was determined
thermal_displacement_matrix_cif – 2D numpy array including the thermal_displacement matrix Ucif format 1st dimension atom types, then compressed thermal displacement matrix will follow U11, U22, U33, U23, U13, U12 (xx, yy, zz, yz, xz, xy) convention similar to “thermal_displacement_matrices.yaml” in phonopy.
- property B: ndarray[source]
Computation as described in R. W. Grosse-Kunstleve, P. D. Adams, J Appl Cryst 2002, 35, 477-480.
- Returns:
First dimension are the atoms in the structure.
- Return type:
np.array
- property U1U2U3: list[source]
Computation as described in R. W. Grosse-Kunstleve, P. D. Adams, J Appl Cryst 2002, 35, 477-480.
- Returns:
eigenvalues of Ucart. First dimension are the atoms in the structure.
- Return type:
np.array
- property Ucif: ndarray[source]
Computation as described in R. W. Grosse-Kunstleve, P. D. Adams, J Appl Cryst 2002, 35, 477-480.
- Returns:
Ucif as array. First dimension are the atoms in the structure.
- Return type:
np.array
- property Ustar: ndarray[source]
Computation as described in R. W. Grosse-Kunstleve, P. D. Adams, J Appl Cryst 2002, 35, 477-480.
- Returns:
Ustar as array. First dimension are the atoms in the structure.
- Return type:
np.array
- property beta: list[source]
Computation as described in R. W. Grosse-Kunstleve, P. D. Adams, J Appl Cryst 2002, 35, 477-480.
- Returns:
First dimension are the atoms in the structure.
- Return type:
np.array
- compute_directionality_quality_criterion(other: ThermalDisplacementMatrices) list[dict[str, ArrayLike]] [source]
Will compute directionality of prolate displacement ellipsoids as described in https://doi.org/10.1039/C9CE00794F with the earlier implementation: https://github.com/damMroz/Angle/.
- Parameters:
other – ThermalDisplacementMatrix
compared (please make sure that the order of the atoms in both objects that are)
Otherwise (is the same.)
results (this analysis will deliver wrong)
- Returns:
will return a list including dicts for each atom that include “vector0” (largest principal axes of self object),
- ”vector1” (largest principal axes of the other object), “angle” between both axes,
These vectors can then, for example, be drawn into the structure with VESTA. Vectors are given in Cartesian coordinates
- classmethod from_Ucif(thermal_displacement_matrix_cif: ArrayLike[ArrayLike], structure: Structure, temperature: float | None = None) Self [source]
Starting from a numpy array, it will convert Ucif values into Ucart values and initialize the class.
- Parameters:
thermal_displacement_matrix_cif – np.array, first dimension are the atoms, then reduced form of thermal displacement matrix will follow Order as above: U11, U22, U33, U23, U13, U12
structure – Structure object
temperature – float Corresponding temperature
- Returns:
ThermalDisplacementMatrices
- static from_cif_P1(filename: str) list[ThermalDisplacementMatrices] [source]
Reads a CIF with P1 symmetry including positions and ADPs. Currently, no check of symmetry is performed as CifParser methods cannot be easily reused.
- Parameters:
filename – Filename of the CIF.
- Returns:
ThermalDisplacementMatrices
- classmethod from_structure_with_site_properties_Ucif(structure: Structure, temperature: float | None = None) Self [source]
Will create this object with the help of a structure with site properties.
- Parameters:
structure – Structure object including U11_cif, U22_cif, U33_cif, U23_cif, U13_cif, U12_cif as site
properties
temperature – temperature for Ucif data
- Returns:
ThermalDisplacementMatrices
- static get_full_matrix(thermal_displacement: ArrayLike[ArrayLike]) np.ndarray[np.ndarray] [source]
Transfers the reduced matrix to the full matrix (order of reduced matrix U11, U22, U33, U23, U13, U12).
- Parameters:
thermal_displacement – 2d numpy array, first dimension are the atoms
- Returns:
3d numpy array including thermal displacements, first dimensions are the atoms
- static get_reduced_matrix(thermal_displacement: ArrayLike[ArrayLike]) np.ndarray[np.ndarray] [source]
Transfers the full matrix to reduced matrix (order of reduced matrix U11, U22, U33, U23, U13, U12).
- Parameters:
thermal_displacement – 2d numpy array, first dimension are the atoms
- Returns:
3d numpy array including thermal displacements, first dimensions are the atoms
- property ratio_prolate: ndarray[source]
This will compute ratio between largest and smallest eigenvalue of Ucart.
- to_structure_with_site_properties_Ucif() Structure [source]
Transfers this object into a structure with site properties (Ucif). This is useful for sorting the atoms in the structure including site properties. e.g. with code like this: def sort_order(site):
return [site.specie.X, site.frac_coords[0], site.frac_coords[1], site.frac_coords[2]]
new_structure0 = Structure.from_sites(sorted(structure0, key=sort_order)).
- Returns:
Structure
- visualize_directionality_quality_criterion(other: ThermalDisplacementMatrices, filename: str | PathLike = 'visualization.vesta', which_structure: Literal[0, 1] = 0) None [source]
Will create a VESTA file for visualization of the directionality criterion.
- Parameters:
other – ThermalDisplacementMatrices
filename – Filename of the VESTA file
which_structure – 0 means structure of the self object will be used, 1 means structure of the other object will be used