API reference

Contents

API reference#

Auto-generated signatures and docstrings for the public phonopy Python API. For a tutorial-style guide, see Phonopy API for Python.

Top-level functions#

phonopy.load

Create a Phonopy instance from parameters and/or input files.

phonopy.load(phonopy_yaml=None, supercell_matrix=None, primitive_matrix='auto', is_nac=True, calculator=None, unitcell=None, supercell=None, nac_params=None, unitcell_filename=None, supercell_filename=None, born_filename=None, force_sets_filename=None, force_constants_filename=None, fc_calculator=None, fc_calculator_options=None, factor=None, produce_fc=True, is_symmetry=True, symmetrize_fc=True, is_compact_fc=True, use_pypolymlp=False, mlp_params=None, use_SNF_supercell=False, symprec=1e-05, log_level=0, lang='Rust')[source]#

Create a Phonopy instance from parameters and/or input files.

A “phonopy_yaml”-like file is parsed unless crystal structure information is given through unitcell_filename, supercell_filename, unitcell (PhonopyAtoms-like), or supercell (PhonopyAtoms-like). Even when a “phonopy_yaml”-like file is parsed, parameters other than crystal structure can be overwritten.

Phonopy’s default FORCE_SETS and BORN files are parsed when they are found in the current directory and the corresponding data are not yet provided by other means.

Crystal structure

Means to provide crystal structure(s), in order of priority:

  1. unitcell_filename (with supercell_matrix)

  2. supercell_filename

  3. unitcell (with supercell_matrix)

  4. supercell

  5. phonopy_yaml

Force sets or force constants

Optional. Means to provide information used to generate force constants, in order of priority:

  1. force_constants_filename

  2. force_sets_filename

  3. phonopy_yaml if force constants are found in it.

  4. phonopy_yaml if forces are found in phonopy_yaml.dataset.

  5. FORCE_CONSTANTS is searched in the current directory.

  6. force_constants.hdf5 is searched in the current directory.

  7. FORCE_SETS is searched in the current directory.

When both 3 and 4 are satisfied but not the others, force constants and dataset are both stored on the Phonopy instance, but force constants are not produced from the dataset.

Parameters for non-analytical term correction (NAC)

Optional. Means to provide NAC parameters, in order of priority:

  1. born_filename

  2. nac_params

  3. phonopy_yaml.nac_params if present and is_nac=True.

  4. BORN is searched in the current directory when is_nac=True.

Parameters:
  • phonopy_yaml (str, os.PathLike, IO, optional) – Filename of a “phonopy.yaml”-like file (str / bytes) or a file-pointer-like object. If given, the data in the file are parsed. Default is None.

  • supercell_matrix (array_like, optional) – Supercell matrix multiplied to the input cell basis vectors. shape=(3,) or (3, 3); the former is interpreted as a diagonal matrix. dtype=int. Default is the identity matrix.

  • primitive_matrix (array_like or str, optional) – Primitive matrix multiplied to the input cell basis vectors. Default is 'auto', which guesses the primitive matrix from crystal symmetry. None is treated the same as 'auto'. To use the unit cell as the primitive cell (identity transformation), pass 'P'. For array_like input, shape=(3, 3), dtype=float. When 'F', 'I', 'A', 'C', or 'R' is given instead of a 3x3 matrix, the primitive matrix for the character found at https://spglib.github.io/spglib/definition.html is used. When a “phonopy.yaml”-like file is loaded and it carries a primitive_matrix, that value takes priority over the default 'auto'.

  • is_nac (bool, optional) – If True, look for a BORN file. If False, NAC is turned off. Default is True.

  • calculator (str, optional) – Calculator used for computing forces. This is used to switch the set of physical units. Default is None, which is equivalent to "vasp".

  • unitcell (PhonopyAtoms, optional) – Input unit cell. Default is None.

  • supercell (PhonopyAtoms, optional) – Input supercell. When given, the default value of primitive_matrix is set to 'auto' (can be overwritten), and supercell_matrix is ignored. Default is None.

  • nac_params (dict, optional) –

    Parameters required for non-analytical term correction. Default is None. Expected structure:

    {'born': Born effective charges
             (array_like, shape=(primitive cell atoms, 3, 3),
              dtype=float),
     'dielectric': Dielectric constant matrix
                   (array_like, shape=(3, 3), dtype=float),
     'factor': unit conversion factor (float)}
    

  • unitcell_filename (str, optional) – Input unit cell filename. Default is None.

  • supercell_filename (str, optional) – Input supercell filename. When specified, supercell_matrix is ignored. Default is None.

  • born_filename (str, optional) – Filename of a BORN-format file containing non-analytical term correction parameters.

  • force_sets_filename (str, optional) – Filename of a file in FORCE_SETS format containing sets of forces and displacements. Default is None.

  • force_constants_filename (str, optional) – Filename of a file in FORCE_CONSTANTS or force_constants.hdf5 format containing force constants. Default is None.

  • fc_calculator ({"traditional", "symfc", "alm", None}, optional) – Force constants calculator backend. Default is None.

  • fc_calculator_options (str, optional) – Optional parameters passed to the external fc-calculator as a single text string. Parsing rules depend on the calculator. For alm, each parameter is split by ',' and each key-value pair is written as 'key = value'.

  • factor (float, optional) – Deprecated. The conversion factor is selected based on calculator.

  • produce_fc (bool, optional) – When False, force constants are not calculated from the dataset of displacements and forces even if the dataset exists. Default is True.

  • is_symmetry (bool, optional) – When False, crystal symmetry (except for lattice translation) is not considered. Default is True.

  • symmetrize_fc (bool, optional) – When False, force constants are not symmetrized when creating them from displacements and forces. Default is True.

  • is_compact_fc (bool, optional) – When True, force constants are produced with shape=(primitive, supercell, 3, 3); when False, with shape=(supercell, supercell, 3, 3). supercell and primitive indicate the number of atoms in those cells. Default is True.

  • use_pypolymlp (bool, optional) – Use pypolymlp for generating force constants. Default is False.

  • mlp_params (dict, optional) – A set of parameters used by machine-learning potentials.

  • use_SNF_supercell (bool, optional) – Build the supercell with the SNF algorithm when True. Default is False. The SNF algorithm is faster than the original one, but the order of atoms in the supercell can be different. Backward compatibility with old data (e.g., force constants) is therefore not guaranteed.

  • symprec (float, optional) – Tolerance used to find crystal symmetry. Default is 1e-5.

  • log_level (int, optional) – Verbosity control. Default is 0.

  • lang (Literal["C", "Rust"], optional) – Backend implementation for compute-heavy kernels. "C" uses the existing C extension; "Rust" selects the experimental phonors backend. Default is "Rust".

Return type:

Phonopy

Phonopy class#

class phonopy.Phonopy(unitcell, supercell_matrix=None, primitive_matrix='auto', factor=None, group_velocity_delta_q=None, symprec=1e-05, is_symmetry=True, use_SNF_supercell=False, hermitianize_dynamical_matrix=True, calculator=None, log_level=0, lang='Rust')[source]#

Bases: object

Phonopy main API.

A Phonopy instance is created from a unit cell and a supercell matrix. It manages displacement generation, force-constant construction, and derived phonon quantities such as band structure, mesh sampling, DOS, thermal properties, group velocity, irreducible representations, modulations, dynamic structure factor, and finite- temperature random displacements.

Most attributes are exposed as @property accessors documented individually below. See Phonopy API for Python for a tutorial-style overview of the typical workflow.

Examples

>>> import numpy as np
>>> from phonopy import Phonopy
>>> from phonopy.structure.atoms import PhonopyAtoms
>>> a = 5.404
>>> unitcell = PhonopyAtoms(
...     symbols=["Si"] * 8,
...     cell=np.eye(3) * a,
...     scaled_positions=[
...         [0, 0, 0], [0, 0.5, 0.5], [0.5, 0, 0.5], [0.5, 0.5, 0],
...         [0.25, 0.25, 0.25], [0.25, 0.75, 0.75],
...         [0.75, 0.25, 0.75], [0.75, 0.75, 0.25],
...     ],
... )
>>> phonon = Phonopy(unitcell, supercell_matrix=[2, 2, 2])
>>> phonon.generate_displacements(distance=0.03)
>>> # Obtain forces by running an external calculator on
>>> # phonon.supercells_with_displacements, then:
>>> # phonon.forces = sets_of_forces
>>> # phonon.produce_force_constants()
>>> # phonon.run_mesh([20, 20, 20])
>>> # phonon.run_thermal_properties(t_step=10, t_max=1000, t_min=0)
Parameters:
  • unitcell (PhonopyAtoms)

  • supercell_matrix (Sequence[int] | Sequence[Sequence[int]] | NDArray | None)

  • primitive_matrix (Literal['P', 'F', 'I', 'A', 'C', 'R', 'auto'] | Sequence[Sequence[float]] | NDArray | None)

  • factor (float | None)

  • group_velocity_delta_q (float | None)

  • symprec (float)

  • is_symmetry (bool)

  • use_SNF_supercell (bool)

  • hermitianize_dynamical_matrix (bool)

  • calculator (str | None)

  • log_level (int)

  • lang (Literal['C', 'Rust'])

property version: str#

Return phonopy release version number.

property primitive: Primitive#

Return primitive cell.

property unitcell: PhonopyAtoms#

Return input unit cell.

property supercell: Supercell#

Return supercell.

property symmetry: Symmetry#

Return symmetry of the supercell.

property primitive_symmetry: Symmetry#

Return symmetry of the primitive cell.

property supercell_matrix: ndarray[tuple[Any, ...], dtype[int64]]#

Return transformation matrix to supercell from unit cell.

Supercell matrix with respect to the unit cell. shape=(3, 3), dtype='int64', order='C'.

property primitive_matrix: ndarray[tuple[Any, ...], dtype[float64]]#

Return transformation matrix to primitive cell from unit cell.

Primitive matrix with respect to the unit cell. shape=(3, 3), dtype='double', order='C'.

property unit_conversion_factor: float#

Return phonon frequency unit conversion factor.

This factor converts sqrt(<force> / <distance> / <AMU>) / 2pi / 1e12 to another preferred phonon frequency unit. It should convert to THz (ordinary frequency) when calculating various phonon properties that assume input phonon frequencies are in THz units. When only frequencies are necessary as output, this factor may be used to get results in other units. The default frequency conversion factor is the one to THz for displacements in Angstroms and forces in eV/Angstrom.

property calculator: str | None#

Return calculator name such as 'vasp', 'qe', etc.

property lang: Literal['C', 'Rust']#

Return the selected backend implementation.

Literal[“C”, “Rust”]

“C” uses the C extension; “Rust” uses the experimental phonors backend.

property dataset: Type1DisplacementDataset | Type2DisplacementDataset | None#

Return displacement-force dataset.

Dataset containing information of displacements in supercells. This optionally contains energies and forces of respective supercells. The format is either one of two types.

Type 1. One atomic displacement in each supercell:

{'natom': number of atoms in supercell,
 'first_atoms': [
   {'number': atom index of displaced atom,
    'displacement': displacement in Cartesian coordinates,
    'forces': forces on atoms in supercell,
    'supercell_energy': energy of supercell},
   {...}, ...]}

Elements of the list accessed by 'first_atoms' correspond to each displaced supercell. Each displaced supercell contains only one displacement. dict['first_atoms']['forces'] gives atomic forces in each displaced supercell.

Type 2. All atomic displacements in each supercell:

{'displacements': ndarray, dtype='double', order='C',
                  shape=(supercells, natom, 3),
 'forces': ndarray, dtype='double', order='C',
                  shape=(supercells, natom, 3),
 'supercell_energies': ndarray, dtype='double'}

To set in type 2, displacements and forces can be given by numpy arrays with different shape that can be reshaped to (supercells, natom, 3).

property mlp_dataset: Type2DisplacementDataset | None#

Return displacement-force dataset used to train an MLP.

The supercell matrix matches that of the usual displacement-force dataset. Only the type-2 format is supported; the dict must contain "displacements", "forces", and "supercell_energies".

property mlp: PhonopyMLP | None#

Setter and getter of the PhonopyMLP machine-learning potential.

property displacements: ndarray[tuple[Any, ...], dtype[float64]] | list#

Getter and setter of displacements in supercells.

There are two types of displacement dataset; see the docstring of dataset for the type-1 and type-2 formats. The returned displacements depend on the dataset type:

Type 1 (list of list)

Each inner list has 4 elements, e.g. [32, 0.01, 0.0, 0.0]. The first element is the supercell atom index starting with 0; the remaining three elements give the displacement in Cartesian coordinates.

Type 2 (ndarray)

Displacements of all atoms in all supercells in Cartesian coordinates. shape=(supercells, natom, 3), dtype='double'.

The setter accepts only the type-2 format (shape=(supercells, natom, 3), dtype='double', order='C').

property force_constants: ndarray[tuple[Any, ...], dtype[float64]] | None#

Getter and setter of supercell force constants.

Force constants matrix.

Getter returns an ndarray with one of two shapes:

  • full: shape=(atoms in supercell, atoms in supercell, 3, 3)

  • compact: shape=(atoms in primitive cell, atoms in supercell, 3, 3)

with dtype='double', order='C'.

Setter accepts any array-like. If given as an own-contiguous ndarray with order='C' and dtype='double', an internal copy of the data is avoided and some computational resources are saved. Expected shape is (atoms in supercell, atoms in supercell, 3, 3), dtype='double'.

set_force_constants_zero_with_radius(cutoff_radius)[source]#

Set zero to force constants within cutoff radius.

Parameters:

cutoff_radius (float)

Return type:

None

property supercell_energies: ndarray[tuple[Any, ...], dtype[float64]]#

Return energies of supercells.

Returns:

shape=(len(supercells),), dtype='double'.

Return type:

ndarray

property forces: ndarray[tuple[Any, ...], dtype[float64]]#

Return forces of supercells.

A set of atomic forces in displaced supercells. The order of displaced supercells has to match with that in the displacement dataset.

The getter returns an ndarray and the setter accepts any array-like with:

shape=(supercells with displacements, atoms in supercell, 3),
dtype='double', order='C'.

That is:

[[[f_1x, f_1y, f_1z], [f_2x, f_2y, f_2z], ...],  # first supercell
 [[f_1x, f_1y, f_1z], [f_2x, f_2y, f_2z], ...],  # second supercell
 ...]
property dynamical_matrix: DynamicalMatrix | None#

Return the DynamicalMatrix instance.

This is the dynamical-matrix builder object, not the matrix itself. Call dm.run(q) and then access dm.dynamical_matrix to obtain the matrix at a given q-point.

property nac_params: dict | None#

Getter and setter of parameters for non-analytical term correction.

A dict with the following entries:

'born'ndarray

Born effective charges. shape=(primitive cell atoms, 3, 3), dtype='double', order='C'.

'dielectric'ndarray

Dielectric constant tensor. shape=(3, 3), dtype='double', order='C'.

'factor'float, optional

Unit conversion factor.

'method'str, optional

Method to calculate NAC.

property supercells_with_displacements: list[PhonopyAtoms] | None#

Return supercells with displacements as a list of PhonopyAtoms.

Generated by generate_displacements().

property mesh_numbers: ndarray[tuple[Any, ...], dtype[int64]] | None#

Return sampling mesh numbers in reciprocal space.

shape=(3,), dtype='int64'. None if run_mesh / init_mesh has not been called.

property qpoints: QpointsPhonon | None#

Return QpointsPhonon instance.

property band_structure: BandStructure | None#

Return BandStructure instance.

property group_velocity: GroupVelocity | None#

Return GroupVelocity instance.

property mesh: Mesh | IterMesh | None#

Return Mesh or IterMesh instance.

property random_displacements: RandomDisplacements | None#

Return RandomDisplacements instance.

property dynamic_structure_factor: DynamicStructureFactor | None#

Return DynamicStructureFactor instance.

property thermal_properties: ThermalProperties | None#

Return ThermalProperties instance.

property thermal_displacements: ThermalDisplacements | None#

Return ThermalDisplacements instance.

property thermal_displacement_matrices: ThermalDisplacementMatrices | None#

Return ThermalDisplacementMatrices instance.

property irreps: IrReps | None#

Return IrReps instance.

property moment: PhononMoment | None#

Return PhononMoment instance.

property total_dos: TotalDos | None#

Return TotalDos instance.

property projected_dos: ProjectedDos | None#

Return ProjectedDos instance.

property masses: ndarray[tuple[Any, ...], dtype[float64]]#

Getter and setter of masses of primitive cell atoms.

By setter, masses of supercell and unit cell atoms are also updated.

generate_displacements(distance=None, is_plusminus='auto', is_diagonal=True, is_trigonal=False, number_of_snapshots=None, random_seed=None, temperature=None, cutoff_frequency=None, max_distance=None, number_estimation_factor=None)[source]#

Generate displacement dataset.

There are two modes, finite difference method with systematic displacements and fitting approach between arbitrary displacements and their forces. The default approach is the finite difference method that is built-in phonopy. The fitting approach requires external force constant calculator.

The random displacement supercells are created by setting positive integer values ‘number_of_snapshots’ keyword argument. Unless this is specified, systematic displacements are created for the finite difference method as the default behaviour.

Parameters:
  • distance (float, optional) – Displacement distance. Unit is the same as that used for the crystal structure. Default is 0.01. For random-direction random-distance displacement generation, this value is also used as min_distance: any generated random distance smaller than this value is replaced by this value.

  • is_plusminus ('auto', True, or False, optional) – For each atom, generate displacements in one direction (False), in both directions (True), or in both directions only when symmetry requires it ('auto'). Default is 'auto'.

  • is_diagonal (bool, optional) – When False, displacements are made only along basis vectors; when True, displacements may be off the basis vectors if doing so reduces the number of displacements by symmetry. Default is True.

  • is_trigonal (bool, optional) – Exists only for testing purposes. Default is False.

  • number_of_snapshots (int, "auto", or None, optional) – Number of snapshots of supercells with random displacements. Random displacements are generated by shifting all atoms in random directions by a fixed distance specified by the distance parameter. In other words, all atoms in the supercell are displaced by the same distance in direct space. When "auto", the minimum required number of snapshots is estimated using symfc and then doubled. The default is None.

  • random_seed (int or None, optional) – Random seed for random displacements generation. Default is None.

  • temperature (float or None, optional) – When given, random displacements at this temperature are generated by sampling the canonical ensemble of harmonic oscillators (harmonic phonons). Default is None.

  • cutoff_frequency (float or None, optional) – In random displacements generation from the canonical ensemble of harmonic phonons, the phonon occupation number is used to determine the deviation of the distribution function. To avoid too large deviation, this value is used to exclude phonon modes whose absolute frequencies are smaller than this value. Default is None.

  • max_distance (float or None, optional) – In random displacements generation from canonical ensemble of harmonic phonons, displacements larger than max distance are renormalized to the max distance, i.e., a displacement d is shortened by d -> d / norm(d) * max_distance if norm(d) > max_distance. In random displacements without using temperature, this value serves as the maximum distance, while the distance parameter sets the minimum distance. The displacement distance is randomly sampled from a uniform distribution between these two bounds.

  • number_estimation_factor (float, optional) – This factor multiplies the number of snapshots estimated by symfc when number_of_snapshots is set to “auto”. Default is None, which sets this factor to 8 when max_distance is specified, otherwise 4.

Return type:

None

produce_force_constants(forces=None, calculate_full_force_constants=True, fc_calculator=None, fc_calculator_options=None, show_drift=True, fc_calculator_log_level=None)[source]#

Compute supercell force constants from forces-displacements dataset.

Supercell force constants are computed from forces and displacements. As the default behaviour, those stored in dataset are used. But with setting forces, this set of forces and the set of displacements stored in the dataset are used for the computation.

Parameters:
  • forces (array_like, optional) – Deprecated. Use the forces setter instead. Default is None.

  • calculate_full_force_constants (bool, optional) – When True, the full force-constants matrix is stored. When False, the compact force-constants matrix is stored. See the docstring of force_constants for details. Default is True.

  • fc_calculator ({"traditional", "symfc", "alm", None}, optional) – Force constants calculator backend. "traditional" uses phonopy’s built-in least-squares fit. "symfc" and "alm" delegate to external packages. Default is None (use the traditional backend).

  • fc_calculator_options (str, optional) – Backend-specific options string passed to the chosen fc_calculator. See the docstring of phonopy.interface.fc_calculator.get_fc2(). Default is None.

  • show_drift (bool, optional) – Display residual translational drift of force constants after computation. Default is True.

  • fc_calculator_log_level (int, optional) – Log level for the force-constants calculator. Default is None (use the Phonopy instance’s log_level).

Return type:

None

symmetrize_force_constants(level=1, show_drift=True, use_symfc_projector=False)[source]#

Symmetrize force constants.

Two schemes are available.

  • Default (use_symfc_projector=False): translational and permutation symmetries are applied successively, not simultaneously. The resulting force constants can break space-group symmetry slightly.

  • use_symfc_projector=True: the symfc projector imposes space-group, translational, and permutation symmetries simultaneously in a single shot.

Parameters:
  • level (int, optional) – Number of times the successive (translation -> permutation) application is repeated. Only used when use_symfc_projector=False. Default is 1.

  • show_drift (bool, optional) – Display residual drift when True. Default is True.

  • use_symfc_projector (bool, optional) – If True, force constants are symmetrized by the symfc projector instead of the traditional approach. Default is False.

Return type:

None

symmetrize_force_constants_by_space_group(show_drift=True)[source]#

Symmetrize force constants using space group operations.

Space group operations except for pure translations are applied to force constants.

Parameters:

show_drift (bool, optional) – Drift forces are displayed when True. Default is True.

Return type:

None

develop_mlp(params=None, test_size=0.1, log_level=None)[source]#

Develop machine learning potential.

Parameters:
  • params (PypolymlpParams or dict, optional) – Parameters for developing MLP. Default is None. When dict is given, PypolymlpParams instance is created from the dict.

  • test_size (float, optional) – Training and test data are split by this ratio. test_size=0.1 means the first 90% of the data is used for training and the rest is used for test. Default is 0.1.

  • log_level (int | None)

Return type:

None

save_mlp(filename=None)[source]#

Save machine learning potential.

Parameters:

filename (str | PathLike | None)

Return type:

None

load_mlp(filename=None)[source]#

Load machine learning potential.

Parameters:

filename (str | PathLike | None)

Return type:

None

evaluate_mlp()[source]#

Evaluate machine learning potential.

This method calculates the supercell energies and forces from the MLP for the displacements in self._dataset of type 2. The results are stored in self._dataset.

The displacements may be generated by the produce_force_constants method with number_of_snapshots > 0. With MLP, a small distance parameter, such as 0.01, can be numerically stable for the computation of force constants.

get_dynamical_matrix_at_q(q)[source]#

Calculate dynamical matrix at a given q-point.

Parameters:

q (array_like) – A q-vector. shape=(3,), dtype='double'.

Returns:

Dynamical matrix. shape=(bands, bands), complex dtype ("c%d" % (np.dtype('double').itemsize * 2)), order='C'.

Return type:

ndarray

get_frequencies(q)[source]#

Calculate phonon frequencies at a given q-point.

Parameters:

q (array_like) – A q-vector. shape=(3,), dtype='double'.

Returns:

Phonon frequencies. Imaginary frequencies are represented by negative real numbers. shape=(bands,), dtype='double'.

Return type:

ndarray

get_frequencies_with_eigenvectors(q)[source]#

Calculate phonon frequencies and eigenvectors at a given q-point.

Parameters:

q (array_like) – A q-vector. shape=(3,).

Returns:

  • frequencies (ndarray) – Phonon frequencies. Imaginary frequencies are represented by negative real numbers. shape=(bands,), dtype='double', order='C'.

  • eigenvectors (ndarray) – Phonon eigenvectors. shape=(bands, bands), complex dtype ("c%d" % (np.dtype('double').itemsize * 2)), order='C'.

Return type:

tuple[ndarray[tuple[Any, …], dtype[float64]], ndarray[tuple[Any, …], dtype[complex128]]]

run_band_structure(paths, with_eigenvectors=False, with_group_velocities=False, is_band_connection=False, path_connections=None, labels=None, is_legacy_plot=False)[source]#

Run phonon band structure calculation.

Parameters:
  • paths (list of array_like) – Sets of q-points defining each band path. The number of q-points can differ between paths. Each array has shape (qpoints, 3).

  • with_eigenvectors (bool, optional) – Whether eigenvectors are calculated. Default is False.

  • with_group_velocities (bool, optional) – Whether group velocities are calculated. Default is False.

  • is_band_connection (bool, optional) – Whether to connect bands across neighboring q-points by comparing the similarity of their eigenvectors. This sometimes fails. Default is False.

  • path_connections (list of bool, optional) – Used only when plotting; indicates whether each path is connected to the next path (i.e., False means there is a jump of q-points between them). The number of elements matches that of paths. Default is None.

  • labels (list of str, optional) – Used only when plotting; labels of the end points of each path. The number of labels equals (2 - np.array(path_connections)).sum().

  • is_legacy_plot (bool, optional) – Use the old-style band-structure plot. Default is False.

Return type:

None

get_band_structure_dict()[source]#

Return calculated band structures.

Returns:

Keys are qpoints, distances, frequencies, eigenvectors, and group_velocities. Each value is a list containing the property along one band path. The number of q-points along one path can be different from that of other paths. Each per-path entry is an ndarray:

qpoints[i]ndarray

q-points in reduced coordinates of reciprocal space without 2 pi. shape=(q-points, 3), dtype='double'.

distances[i]ndarray

Distances in reciprocal space along paths. shape=(q-points,), dtype='double'.

frequencies[i]ndarray

Phonon frequencies. Imaginary frequencies are represented by negative real numbers. shape=(q-points, bands), dtype='double'.

eigenvectors[i]ndarray

Phonon eigenvectors. None if eigenvectors are not stored. shape=(q-points, bands, bands), dtype=complex ("c%d" % (np.dtype('double').itemsize * 2)), order='C'.

group_velocities[i]ndarray

Phonon group velocities. None if group velocities are not calculated. shape=(q-points, bands, 3), dtype='double'.

Return type:

dict

auto_band_structure(npoints=101, with_eigenvectors=False, with_group_velocities=False, plot=False, write_yaml=False, filename='band.yaml')[source]#

Conveniently calculate and draw band structure.

See the docstring of run_band_structure() for the parameters with_eigenvectors (default False) and with_group_velocities (default False).

Parameters:
  • npoints (int, optional) – Number of q-points in each segment of the band-structure paths. The number includes end points. Default is 101.

  • plot (bool, optional) – When True, band structure is plotted using matplotlib and the matplotlib module (plt) is returned. To watch the result, usually show() has to be called. Default is False.

  • write_yaml (bool, optional) – When True, a band.yaml like file is written out. The file name can be specified with the filename parameter. Default is False.

  • filename (str, optional) – File name used to write the band.yaml like file. Default is band.yaml.

  • with_eigenvectors (bool)

  • with_group_velocities (bool)

Return type:

Any | None

plot_band_structure()[source]#

Plot calculated band structure.

Returns:

The matplotlib.pyplot module. Call .show() on it to display the figure.

Return type:

matplotlib.pyplot

write_hdf5_band_structure(comment=None, filename='band.hdf5', compression=None)[source]#

Write band structure in hdf5 format.

Parameters:
  • comment (dict, optional) – Items are stored in hdf5 file in the way of key-value pair.

  • filename (str, optional) – Default is band.hdf5.

  • compression (Literal['gzip', 'lzf'] | int | None)

Return type:

None

write_yaml_band_structure(comment=None, filename=None, compression=None)[source]#

Write band structure in yaml.

Parameters:
  • comment (dict) – Data structure dumped in YAML and the dumped YAML text is put at the beginning of the file.

  • filename (str) – Default filename is ‘band.yaml’ when compression=None. With compression, an extension of filename is added such as ‘band.yaml.xz’.

  • compression (None, 'gzip', or 'lzma') – None gives a plain text file. 'gzip' and 'lzma' compress the yaml text with the respective compression method.

Return type:

None

init_mesh(mesh=100.0, shift=None, is_time_reversal=True, is_mesh_symmetry=True, with_eigenvectors=False, with_group_velocities=False, is_gamma_center=False, use_iter_mesh=False)[source]#

Initialize mesh sampling phonon calculation without starting to run.

Phonon calculation starts explicitly with calling Mesh.run() or implicitly with accessing getters of Mesh instance, e.g., Mesh.frequencies.

Parameters:
  • mesh (array_like or float, optional) – Mesh numbers along a, b, c axes when array_like object is given. dtype='int64', shape=(3,). When a float value is given, a uniform mesh is generated following the VASP convention by N = max(1, nint(l * norm(a*))), where nint is the function that returns the nearest integer and a* is each reciprocal basis vector. In this case, is_gamma_center=True is enforced. Default value is 100.0.

  • shift (array_like, optional) – Mesh shifts along a*, b*, c* axes with respect to neighboring grid points from the original mesh (Monkhorst-Pack or Gamma center). 0.5 gives a half-grid shift. Normally 0 or 0.5 is given; otherwise q-point symmetry search is not performed. Default is None (no additional shift). shape=(3,), dtype='double'.

  • is_time_reversal (bool, optional) – Whether to include time-reversal symmetry in the symmetry search. Default is True.

  • is_mesh_symmetry (bool, optional) – Whether mesh symmetry search is performed. Default is True.

  • with_eigenvectors (bool, optional) – Store eigenvectors when True. Default is False.

  • with_group_velocities (bool, optional) – Calculate group velocities when True. Default is False.

  • is_gamma_center (bool, optional) – Generate a uniform mesh centered at Gamma instead of using the Monkhorst-Pack scheme. When mesh is given as a float (length measure), this setting is ignored and is_gamma_center=True is enforced. Default is False.

  • use_iter_mesh (bool, optional) – Use IterMesh instead of Mesh so that phonon properties are not stored on the instance, saving memory. Used with ThermalDisplacements and ThermalDisplacementMatrices. Default is False.

Return type:

None

run_mesh(mesh=100.0, shift=None, is_time_reversal=True, is_mesh_symmetry=True, with_eigenvectors=False, with_group_velocities=False, is_gamma_center=False)[source]#

Run mesh sampling phonon calculation.

See the parameter details in Phonopy.init_mesh.

Parameters:
Return type:

None

get_mesh_dict()[source]#

Return phonon properties calculated by mesh sampling.

Returns:

Keys are qpoints, weights, frequencies, eigenvectors, and group_velocities.

qpointsndarray

q-points in reduced coordinates of the reciprocal lattice. shape=(ir-grid points, 3), dtype='double'.

weightsndarray

Geometric q-point weights. The sum equals the number of grid points. shape=(ir-grid points,), dtype='int64'.

frequenciesndarray

Phonon frequencies at ir-grid points. Imaginary frequencies are represented by negative real numbers. shape=(ir-grid points, bands), dtype='double'.

eigenvectorsndarray

Phonon eigenvectors at ir-grid points. See the data structure of np.linalg.eigh. shape=(ir-grid points, bands, bands), complex dtype ("c%d" % (np.dtype('double').itemsize * 2)), order='C'.

group_velocitiesndarray

Phonon group velocities at ir-grid points. shape=(ir-grid points, bands, 3), dtype='double'.

Return type:

dict

write_hdf5_mesh(compression=None)[source]#

Write mesh calculation results in hdf5 format.

Parameters:

compression (Literal['gzip', 'lzf'] | int | None)

Return type:

None

write_yaml_mesh()[source]#

Write mesh calculation results in yaml format.

Return type:

None

plot_band_structure_and_dos(pdos_indices=None)[source]#

Plot band structure and DOS.

Parameters:

pdos_indices (Sequence[Sequence[int]] | None)

Return type:

Any

run_qpoints(q_points, with_eigenvectors=False, with_group_velocities=False, with_dynamical_matrices=False, nac_q_direction=None)[source]#

Run phonon calculation at specified q-points.

Parameters:
  • q_points (array_like) – q-points in reduced coordinates. dtype=’double’, shape=(q-points, 3)

  • with_eigenvectors (bool, optional) – Eigenvectors are stored by setting True. Default False.

  • with_group_velocities (bool, optional) – Group velocities are calculated by setting True. Default is False.

  • with_dynamical_matrices (bool, optional) – Calculated dynamical matrices are stored by setting True. Default is False.

  • nac_q_direction (array_like, optional) – q-point direction from Gamma-point in fractional coordinates of reciprocal basis vectors. Only the direction is used, i.e., q_direction / norm(q_direction) is computed and used. This parameter is activated only at q=(0, 0, 0). shape=(3,), dtype='double'.

Return type:

None

get_qpoints_dict()[source]#

Return calculated phonon properties at q-points.

Returns:

Keys are frequencies, eigenvectors, group_velocities, and dynamical_matrices.

frequenciesndarray

Phonon frequencies. Imaginary frequencies are represented by negative real numbers. shape=(qpoints, bands), dtype='double'.

eigenvectorsndarray or None

Phonon eigenvectors. None if eigenvectors are not stored. shape=(qpoints, bands, bands), complex dtype ("c%d" % (np.dtype('double').itemsize * 2)), order='C'.

group_velocitiesndarray or None

Phonon group velocities. None if group velocities are not calculated. shape=(qpoints, bands, 3), dtype='double'.

dynamical_matricesndarray

Dynamical matrices at q-points. shape=(qpoints, bands, bands), dtype='double'.

Return type:

dict

write_hdf5_qpoints_phonon(compression=None)[source]#

Write phonon properties calculated at q-points in hdf5 format.

Parameters:

compression (Literal['gzip', 'lzf'] | int | None)

Return type:

None

write_yaml_qpoints_phonon()[source]#

Write phonon properties calculated at q-points in yaml format.

Return type:

None

run_total_dos(sigma=None, freq_min=None, freq_max=None, freq_pitch=None, use_tetrahedron_method=True)[source]#

Run total DOS calculation.

Parameters:
  • sigma (float, optional) – Smearing width for the smearing method. Default is None.

  • freq_min (float, optional) – Minimum and maximum frequencies of the frequency range in which DOS is computed, and the sampling interval (freq_pitch). Defaults are None and they are automatically determined.

  • freq_max (float, optional) – Minimum and maximum frequencies of the frequency range in which DOS is computed, and the sampling interval (freq_pitch). Defaults are None and they are automatically determined.

  • freq_pitch (float, optional) – Minimum and maximum frequencies of the frequency range in which DOS is computed, and the sampling interval (freq_pitch). Defaults are None and they are automatically determined.

  • use_tetrahedron_method (bool, optional) – Use the tetrahedron method when True. When sigma is set, the smearing method is used instead. Default is True.

Return type:

None

auto_total_dos(mesh=100.0, is_time_reversal=True, is_mesh_symmetry=True, is_gamma_center=False, plot=False, xlabel=None, ylabel=None, with_tight_frequency_range=False, write_dat=False, filename='total_dos.dat')[source]#

Conveniently calculate and draw total DOS.

Return type:

Any | None

get_total_dos_dict()[source]#

Return total DOS.

Returns:

  • A dictionary with keys of ‘frequency_points’ and ‘total_dos’.

  • Each value of corresponding key is as follows

  • frequency_points (ndarray) – shape=(frequency_sampling_points, ), dtype=’double’

  • total_dos – shape=(frequency_sampling_points, ), dtype=’double’

Return type:

TotalDosDict

set_Debye_frequency(freq_max_fit=None)[source]#

Calculate Debye frequency on top of total DOS.

Parameters:

freq_max_fit (float | None)

Return type:

None

get_Debye_frequency()[source]#

Return Debye frequency.

Return type:

float | None

plot_total_dos(xlabel=None, ylabel=None, with_tight_frequency_range=False)[source]#

Plot total DOS.

Parameters:
  • xlabel (str, optional) – x-label of the plot. Default is None, which puts a default x-label.

  • ylabel (str, optional) – y-label of the plot. Default is None, which puts a default y-label.

  • with_tight_frequency_range (bool, optional) – Plot with a tight frequency range. Default is False.

Return type:

Any

write_total_dos(filename='total_dos.dat')[source]#

Write total DOS to text file.

Parameters:

filename (str | PathLike)

Return type:

None

run_projected_dos(sigma=None, freq_min=None, freq_max=None, freq_pitch=None, use_tetrahedron_method=True, direction=None, xyz_projection=False)[source]#

Run projected DOS calculation.

Parameters:
  • sigma (float, optional) – Smearing width for the smearing method. Default is None.

  • freq_min (float, optional) – Minimum and maximum frequencies of the frequency range in which DOS is computed, and the sampling interval (freq_pitch). Defaults are None and they are automatically determined.

  • freq_max (float, optional) – Minimum and maximum frequencies of the frequency range in which DOS is computed, and the sampling interval (freq_pitch). Defaults are None and they are automatically determined.

  • freq_pitch (float, optional) – Minimum and maximum frequencies of the frequency range in which DOS is computed, and the sampling interval (freq_pitch). Defaults are None and they are automatically determined.

  • use_tetrahedron_method (bool, optional) – Use the tetrahedron method when True. When sigma is set, the smearing method is used instead. Default is True.

  • direction (array_like, optional) – Projection direction given as three values along the primitive cell basis vectors. Default is None (no projection).

  • xyz_projection (bool, optional) – Whether to project along Cartesian directions. Default is False.

Return type:

None

auto_projected_dos(mesh=100.0, is_time_reversal=True, is_gamma_center=False, plot=False, pdos_indices=None, legend=None, legend_prop=None, legend_frameon=True, xlabel=None, ylabel=None, with_tight_frequency_range=False, write_dat=False, filename='projected_dos.dat')[source]#

Conveniently calculate and draw projected DOS.

See the docstring of Phonopy.init_mesh for the parameters mesh (default 100.0), is_time_reversal (default True), and is_gamma_center (default False). See the docstring of Phonopy.plot_projected_dos for pdos_indices, legend, xlabel, ylabel, and with_tight_frequency_range.

Parameters:
  • plot (bool, optional) – With setting True, PDOS is plotted using matplotlib and the matplotlib module (plt) is returned. To watch the result, usually show() has to be called. Default is False.

  • write_dat (bool, optional) – With setting True, a projected_dos.dat like file is written out. The file name can be specified with the filename parameter. Default is False.

  • filename (str, optional) – File name used to write the projected_dos.dat like file. Default is projected_dos.dat.

  • mesh (float | Sequence[int] | ndarray[tuple[Any, ...], dtype[int64]])

  • is_time_reversal (bool)

  • is_gamma_center (bool)

  • pdos_indices (Sequence[Sequence[int]] | None)

  • legend (Sequence[str] | None)

  • legend_prop (dict | None)

  • legend_frameon (bool)

  • xlabel (str | None)

  • ylabel (str | None)

  • with_tight_frequency_range (bool)

Return type:

Any | None

get_projected_dos_dict()[source]#

Return projected DOS.

Projection is done to atoms and may be also done along directions depending on the parameters at run_projected_dos.

Returns:

  • A dictionary with keys of ‘frequency_points’ and ‘projected_dos’.

  • Each value of corresponding key is as follows

  • frequency_points (ndarray) – shape=(frequency_sampling_points, ), dtype=’double’

  • projected_dos – shape=(projections, frequency_sampling_points), dtype=’double’

Return type:

ProjectedDosDict

plot_projected_dos(pdos_indices=None, legend=None, legend_prop=None, legend_frameon=True, xlabel=None, ylabel=None, with_tight_frequency_range=False)[source]#

Plot projected DOS.

Parameters:
  • pdos_indices (list of list, optional) – Sets of indices of atoms whose projected DOS are summed over. The indices start with 0. An example is pdos_indices=[[0, 1], [2, 3, 4, 5]]. Default is None, which means pdos_indices=[[i] for i in range(natom)].

  • legend (list of instances such as str or int, optional) – The str(instance) are shown in legend. It has to be len(pdos_indices)==len(legend). Default is None. When None, legend is not shown.

  • legend_prop (dict, optional) – Legend properties of matplotlib. Default is None.

  • legend_frameon (bool, optional) – Legend with frame or not. Default is True.

  • xlabel (str, optional) – x-label of plot. Default is None, which puts a default x-label.

  • ylabel (str, optional) – y-label of plot. Default is None, which puts a default y-label.

  • with_tight_frequency_range (bool, optional) – Plot with tight frequency range. Default is False.

Return type:

Any

write_projected_dos(filename='projected_dos.dat')[source]#

Write projected DOS to text file.

Parameters:

filename (str | PathLike)

Return type:

None

run_thermal_properties(t_min=0, t_max=1000, t_step=10, temperatures=None, cutoff_frequency=None, pretend_real=False, band_indices=None, is_projection=False, classical=False)[source]#

Run calculation of thermal properties at constant volume.

In phonopy, imaginary frequencies are represented as negative real value. Under this situation, cutoff_frequency is used to ignore phonon modes that have frequencies less than cutoff_frequency.

Parameters:
  • t_min (float, optional) – Minimum and maximum temperatures and the interval in this temperature range. Default values are 0, 1000, and 10.

  • t_max (float, optional) – Minimum and maximum temperatures and the interval in this temperature range. Default values are 0, 1000, and 10.

  • t_step (float, optional) – Minimum and maximum temperatures and the interval in this temperature range. Default values are 0, 1000, and 10.

  • temperatures (array_like, optional) – Temperature points where thermal properties are calculated. When this is set, t_min, t_max, and t_step are ignored.

  • cutoff_frequency (float, optional) – Ignore phonon modes whose frequencies are smaller than this value. Default is None, which gives cutoff frequency as zero.

  • pretend_real (bool, optional) – Use absolute value of phonon frequency when True. Default is False.

  • band_indices (array_like, optional) – Band indices starting with 0. Normally the numbers correspond to phonon bands in ascending order of phonon frequencies. Thermal properties are calculated only including specified bands. Note that use of this results in unphysical values, and it is not recommended to use this feature. Default is None.

  • is_projection (bool, optional) – When True, fractions of squared eigenvector elements are multiplied with the mode thermal-property quantities at the respective phonon modes. Note that using this results in unphysical values, so this feature is not recommended. Default is False.

  • classical (bool, optional) – If True, use classical statistics; if False, use quantum statistics. Default is False.

Return type:

None

get_thermal_properties_dict()[source]#

Return thermal properties.

Returns:

  • A dictionary of thermal properties with keys of ‘temperatures’,

  • ’free_energy’, ‘entropy’, and ‘heat_capacity’.

  • Each value of corresponding key is as follows

  • temperatures (ndarray) – shape=(temperatures, ), dtype=’double’

  • free_energy (ndarray) – shape=(temperatures, ), dtype=’double’

  • entropy (ndarray) – shape=(temperatures, ), dtype=’double’

  • heat_capacity (ndarray) – shape=(temperatures, ), dtype=’double’

Return type:

ThermalPropertiesDict

plot_thermal_properties(xlabel=None, ylabel=None, with_grid=True, divide_by_Z=False, legend_style='normal')[source]#

Plot thermal properties.

Parameters:
  • xlabel (str, optional) – Label used for x-axis.

  • ylabel (str, optional) – Label used for y-axis.

  • with_grid (bool, optional) – With grid or not. Default is True.

  • divide_by_Z (bool, optional) – Divide thermal properties by number of formula units of primitive cell. Default is False.

  • legend_style (str, optional) – “normal”, “compact”, None. None will not show legend.

Return type:

Any

write_yaml_thermal_properties(filename='thermal_properties.yaml')[source]#

Write thermal properties in yaml format.

Parameters:

filename (str | PathLike)

Return type:

None

run_thermal_displacements(t_min=0, t_max=1000, t_step=10, temperatures=None, direction=None, freq_min=None, freq_max=None)[source]#

Run thermal displacements calculation.

Parameters:
  • t_min (float, optional) – Minimum and maximum temperatures and the interval in this temperature range. Default values are 0, 1000, and 10.

  • t_max (float, optional) – Minimum and maximum temperatures and the interval in this temperature range. Default values are 0, 1000, and 10.

  • t_step (float, optional) – Minimum and maximum temperatures and the interval in this temperature range. Default values are 0, 1000, and 10.

  • temperatures (array_like, optional) – Temperature points where thermal properties are calculated. When this is set, t_min, t_max, and t_step are ignored.

  • direction (array_like, optional) – Projection direction in reduced coordinates. shape=(3,), dtype='double'. Default is None (no projection).

  • freq_min (float, optional) – Only phonon frequencies between freq_min and freq_max are included. Default is None (all phonons).

  • freq_max (float, optional) – Only phonon frequencies between freq_min and freq_max are included. Default is None (all phonons).

Return type:

None

get_thermal_displacements_dict()[source]#

Return thermal displacements.

Return type:

dict

plot_thermal_displacements(is_legend=False)[source]#

Plot thermal displacements.

Parameters:

is_legend (bool)

Return type:

Any

write_yaml_thermal_displacements()[source]#

Write thermal displacements in yaml format.

Return type:

None

run_thermal_displacement_matrices(t_min=0, t_max=1000, t_step=10, temperatures=None, freq_min=None, freq_max=None)[source]#

Run thermal displacement matrices calculation.

Parameters:
  • t_min (float, optional) – Minimum and maximum temperatures and the interval in this temperature range. Default values are 0, 1000, and 10.

  • t_max (float, optional) – Minimum and maximum temperatures and the interval in this temperature range. Default values are 0, 1000, and 10.

  • t_step (float, optional) – Minimum and maximum temperatures and the interval in this temperature range. Default values are 0, 1000, and 10.

  • freq_min (float, optional) – Phonon frequencies larger than freq_min and smaller than freq_max are included. Default is None, i.e., all phonons.

  • freq_max (float, optional) – Phonon frequencies larger than freq_min and smaller than freq_max are included. Default is None, i.e., all phonons.

  • temperatures (array_like, optional) – Temperature points where thermal properties are calculated. When this is set, t_min, t_max, and t_step are ignored. Default is None.

Return type:

None

get_thermal_displacement_matrices_dict()[source]#

Return thermal displacement matrices.

Return type:

ThermalDisplacementMatricesDict

write_yaml_thermal_displacement_matrices()[source]#

Write thermal displacement matrices in yaml format.

Return type:

None

write_thermal_displacement_matrix_to_cif(temperature_index)[source]#

Write thermal displacement matrices at a temperature in cif.

Parameters:

temperature_index (int)

Return type:

None

write_animation(q_point=None, anime_type='v_sim', band_index=None, amplitude=None, num_div=None, shift=None, filename=None)[source]#

Write atomic modulations in animation format.

Returns:

Output filename.

Return type:

str

Parameters:
run_modulations(dimension, phonon_modes, delta_q=None, derivative_order=None, nac_q_direction=None)[source]#

Generate atomic displacements of phonon modes.

The design of this feature, and thus its API, is not very satisfactory. It should be reconsidered someday in the future.

Parameters:
  • dimension (array_like) – Supercell dimension with respect to the primitive cell. shape=(3,), (3, 3), or (9,), dtype='int64'.

  • phonon_modes (list of phonon mode settings) –

    Each element of the outer list specifies one phonon mode:

    [q-point, band index (int), amplitude (float),
     phase (float)]
    

    The first element is a list representing the q-point in reduced coordinates. The remaining elements are the band index (starting with 0), amplitude, and phase factor.

  • nac_q_direction (array_like) – q-point direction from Gamma-point in fractional coordinates of reciprocal basis vectors. Only the direction is used, i.e., q_direction / norm(q_direction) is computed and used. This parameter is activated only at q=(0, 0, 0). shape=(3,), dtype='double'.

  • delta_q (Sequence[float] | ndarray[tuple[Any, ...], dtype[float64]] | None)

  • derivative_order (int | None)

Return type:

None

get_modulated_supercells()[source]#

Return modulated structures as a list of PhonopyAtoms.

Return type:

list[PhonopyAtoms]

get_modulations_and_supercell()[source]#

Return atomic modulations and the perfect supercell.

Returns:

  • modulations (ndarray) – Atomic modulations of the supercell in Cartesian coordinates.

  • supercell (PhonopyAtoms) – The (unmodulated) supercell.

Return type:

tuple[ndarray[tuple[Any, …], dtype[complex128]], PhonopyAtoms]

write_modulations(calculator=None, optional_structure_info=None)[source]#

Write modulated structures to MPOSCAR files.

Parameters:
  • calculator (str | None)

  • optional_structure_info (StructureInfo | None)

Return type:

None

write_yaml_modulations()[source]#

Write atomic modulations in yaml format.

Return type:

None

set_irreps(q, is_little_cogroup=False, nac_q_direction=None, degeneracy_tolerance=None)[source]#

Identify ir-reps of phonon modes.

The design of this API is not very satisfactory and is expected to be redesigned in the next major versions once the use case of the API for the ir-reps feature becomes clearer.

Parameters:
  • nac_q_direction (array_like) – q-point direction from Gamma-point in fractional coordinates of reciprocal basis vectors. Only the direction is used, i.e., q_direction / norm(q_direction) is computed and used. This parameter is activated only at q=(0, 0, 0). shape=(3,), dtype='double'.

  • q (Sequence[float] | ndarray[tuple[Any, ...], dtype[float64]])

  • is_little_cogroup (bool)

  • degeneracy_tolerance (float | None)

Return type:

None

show_irreps(show_irreps=False)[source]#

Show Ir-reps.

Parameters:

show_irreps (bool)

Return type:

None

write_yaml_irreps(show_irreps=False)[source]#

Write Ir-reps in yaml format.

Parameters:

show_irreps (bool)

Return type:

None

get_group_velocity_at_q(q_point)[source]#

Return group velocity at a q-point.

Parameters:

q_point (Sequence[float] | ndarray[tuple[Any, ...], dtype[float64]])

Return type:

ndarray[tuple[Any, …], dtype[float64]]

run_moment(order=1, is_projection=False, freq_min=None, freq_max=None)[source]#

Run moment calculation.

Parameters:
Return type:

None

get_moment()[source]#

Return moment.

Return type:

float | ndarray[tuple[Any, …], dtype[float64]] | None

init_dynamic_structure_factor(Qpoints, T, atomic_form_factor_func=None, scattering_lengths=None, freq_min=None, freq_max=None)[source]#

Initialize dynamic structure factor calculation.

Call DynamicStructureFactor.run() to start the calculation.

Parameters:
  • Qpoints (array_like) – Q-points in any Brillouin zone. shape=(qpoints, 3), dtype='double'.

  • T (float) – Temperature in K.

  • atomic_form_factor_func (callable) –

    Function that returns the atomic form factor (func below):

    f_params = {
        'Na': [3.148690, 2.594987, 4.073989, 6.046925,
               0.767888, 0.070139, 0.995612, 14.1226457,
               0.968249, 0.217037, 0.045300],
        'Cl': [1.061802, 0.144727, 7.139886, 1.171795,
               6.524271, 19.467656, 2.355626, 60.320301,
               35.829404, 0.000436, -34.916604],
    }
    
    def get_func_AFF(f_params):
        def func(symbol, Q):
            return atomic_form_factor_WK1995(Q, f_params[symbol])
        return func
    

  • scattering_lengths (dict) – Coherent scattering lengths averaged over isotopes and spins. Supposed for INS. For example, {'Na': 3.63, 'Cl': 9.5770}.

  • freq_min (float) – Minimum and maximum phonon frequencies to determine whether phonons are included in the calculation.

  • freq_max (float) – Minimum and maximum phonon frequencies to determine whether phonons are included in the calculation.

Return type:

None

run_dynamic_structure_factor(Qpoints, T, atomic_form_factor_func=None, scattering_lengths=None, freq_min=None, freq_max=None)[source]#

Run dynamic structure factor calculation.

See the detail of parameters at Phonopy.init_dynamic_structure_factor().

Parameters:
Return type:

None

get_dynamic_structure_factor()[source]#

Return dynamic structure factors.

Return type:

tuple[ndarray[tuple[Any, …], dtype[float64]], ndarray[tuple[Any, …], dtype[float64]]]

init_random_displacements(dist_func=None, cutoff_frequency=None, max_distance=None)[source]#

Initialize random displacements at finite temperature.

Parameters:
  • dist_func (str or None, optional) – Harmonic oscillator distribution function: either 'quantum' or 'classical'. Default is None, corresponding to 'quantum'.

  • cutoff_frequency (float or None) – Phonon frequency in THz below which phonons are ignored when generating random displacements. Default is None.

  • max_distance (float or None, optional) – In random displacements generation from canonical ensemble of harmonic phonons, displacements larger than max distance are renormalized to the max distance, i.e., a displacement d is shortened by d -> d / norm(d) * max_distance if norm(d) > max_distance.

Return type:

None

get_random_displacements_at_temperature(temperature, number_of_snapshots, is_plusminus=False, random_seed=None)[source]#

Generate random displacements from phonon structure.

See generate_displacements() for related details.

Parameters:
  • temperature (float) – Temperature in K.

  • number_of_snapshots (int) – Number of snapshots with random displacements to create.

  • is_plusminus (bool, optional) – If True, concatenate the displacements with their opposites, doubling the number of snapshots. Default is False.

  • random_seed (32-bit unsigned int or None, optional) – Random seed. Default is None.

Return type:

ndarray

save(filename='phonopy_params.yaml', settings=None, hdf5_settings=None, compression=False)[source]#

Save phonopy parameters into file.

Parameters:
  • filename (str, optional) – File name. Default is "phonopy_params.yaml".

  • settings (dict, optional) –

    Selects which parameters are written out. Only the entries to be changed from the defaults need to be set. The available parameters and their default settings are:

    {'force_sets': True,
     'displacements': True,
     'force_constants': False,
     'born_effective_charge': True,
     'dielectric_constant': True}
    

    The default settings are updated by {'force_constants': True} when dataset is None and force_constants is not None, unless {'force_constants': False} is specified explicitly.

  • hdf5_settings (dict, optional (to be implemented)) –

    Force constants and force sets are stored in an HDF5 file when activated in the dict. The default filename is the filename of the yaml file with .yaml replaced by .hdf5. Keys:

    {'filename': str,
     'force_constants': bool (default=False),
     'force_sets': bool (default=False)}
    

  • compression (bool or str) – If True, the phonopy_params.yaml like file is compressed by xz. When compression == 'xz', the file is compressed by xz. Default is False.

Returns:

File name of the saved phonopy_params.yaml like file (with an .xz suffix if compressed).

Return type:

str

ph2ph(supercell_matrix, with_nac=False)[source]#

Transform force constants in this Phonopy instance to another shape.

Force constants are Fourier-interpolated. This Phonopy instance must already hold force constants. The init parameters of this instance are copied to the returned instance.

For example, if the current supercell_matrix is [2, 2, 2] and the given supercell_matrix is [4, 4, 4], the existing force constants are Fourier-interpolated by sampling at the commensurate points of the latter supercell, and a new Phonopy instance carrying the interpolated force constants is returned.

Parameters:
  • supercell_matrix (array_like) – Specifies the shape of the new force constants.

  • with_nac (bool, optional) – Use non-analytical term correction (NAC) under the Fourier interpolation: dynamical matrices at commensurate points are computed with NAC, then Fourier-transformed back to force constants of the requested supercell_matrix. NAC parameters are not copied to the returned Phonopy instance.

Returns:

Phonopy instance carrying the init parameters of this instance and the transformed force constants for the given supercell_matrix.

Return type:

Phonopy

copy(log_level=None)[source]#

Copy this Phonopy instance carrying only the init parameters.

Notes

The returned instance is constructed with the same init parameters as this one, but internal state such as force constants, NAC parameters, MLP, etc. is not copied.

Returns:

New Phonopy instance.

Return type:

Phonopy

Parameters:

log_level (int | None)

to_phonopy_yaml(configuration=None, settings=None)[source]#

Return PhonopyYaml class instance with this data.

Parameters:
  • configuration (dict | None)

  • settings (dict | None)

Return type:

PhonopyYaml

PhonopyAtoms class#

class phonopy.structure.atoms.PhonopyAtoms(symbols=None, numbers=None, masses=None, magnetic_moments=None, scaled_positions=None, positions=None, cell=None, species_table=None, species_ids=None)[source]#

Bases: object

Class to represent crystal structure.

Originally this class aimed to be compatible with the ASE Atoms class, but the APIs have since diverged.

A cell is described by basis vectors (cell), per-atom positions (positions / scaled_positions), and per-atom species information. Species can be expressed either as plain chemical symbols (symbols) / atomic numbers (numbers) for ordinary cells, or as a deduplicated table plus per-atom indices (species_table / species_ids) which additionally supports mixed-species sites for the Virtual Crystal Approximation. See PhonopyAtoms class for the tutorial-style overview and per-attribute documentation below for details.

Parameters:
  • symbols (Sequence[str] | None)

  • numbers (Sequence[int] | NDArray[np.int64] | None)

  • masses (Sequence[float] | NDArray[np.double] | None)

  • magnetic_moments (Sequence[float] | Sequence[Sequence[float]] | NDArray[np.double] | None)

  • scaled_positions (Sequence[Sequence[float]] | NDArray[np.double] | None)

  • positions (Sequence[Sequence[float]] | NDArray[np.double] | None)

  • cell (Sequence[Sequence[float]] | NDArray[np.double] | None)

  • species_table (Sequence[_Species] | None)

  • species_ids (Sequence[int] | NDArray[np.int64] | None)

property cell: ndarray[tuple[Any, ...], dtype[float64]]#

Setter and getter of basis vectors.

Basis vectors (a, b, c) given as row vectors. shape=(3, 3), dtype='double', order='C'. For getter, a copy is returned.

property positions: ndarray[tuple[Any, ...], dtype[float64]]#

Setter and getter of positions in Cartesian coordinates.

shape=(natom, 3), dtype='double', order='C'.

property scaled_positions: ndarray[tuple[Any, ...], dtype[float64]]#

Setter and getter of scaled positions.

Positions of atoms in fractional (crystallographic) coordinates. shape=(natom, 3), dtype='double', order='C'. For getter, a copy is returned.

property symbols: list[str]#

Chemical symbols per atom.

Chemical symbol with an appended natural number is allowed, e.g., "Cl1". Mixed-species sites carry a composite label formed by concatenating constituent symbols (e.g., "GeSn").

property species_ids: ndarray[tuple[Any, ...], dtype[int64]]#

Per-atom index into species_table.

Atoms that reference the same species_table entry share the same id; ordinary "Cl" and labeled "Cl1" get different ids, and two atoms differing only in mixture composition or weights also get different ids. Suitable as the types argument for spglib when atoms with the same atomic number but different symbol suffix must be distinguished. shape=(natom,), dtype='int64'. For getter, a copy is returned.

property species_table: list[_Species]#

Deduplicated species table indexed by species_ids.

Each entry is either an ordinary single-element species (with atomic_number set) or a weighted mixture (mixture set to a tuple of (symbol, weight) pairs summing to 1.0). See phonopy.structure.atoms.build_species_table_from_mixtures(). For getter, a shallow list copy is returned; _Species entries are frozen dataclasses so the list elements are themselves immutable.

property numbers: ndarray[tuple[Any, ...], dtype[int64]]#

Atomic numbers per atom in 1..118.

shape=(natom,), dtype='int64'. For getter, a new array is returned.

Raises RuntimeError if the cell contains any mixed-species site, since a mixture has no single atomic number. Use species_ids instead when an opaque per-atom species discriminator is needed.

property has_mixtures: bool#

Return True if the cell contains any mixed-species site.

True when any species in the cell is a weighted mixture of constituents (e.g., a Virtual Crystal Approximation site).

property masses: ndarray[tuple[Any, ...], dtype[float64]]#

Setter and getter of atomic masses.

For a mixed-species site, the mass is the weight-averaged sum over its constituents. shape=(natom,), dtype='double'. For getter, a copy is returned.

property magnetic_moments: ndarray[tuple[Any, ...], dtype[float64]] | None#

Setter and getter of magnetic moments.

shape=(natom,) or (natom, 3), dtype='double', order='C'.

For the setter, the former can also be specified as (natom, 1) (recognized as (natom,)) and the latter as (natom * 3,) (converted to (natom, 3)). For getter, a copy is returned.

property volume: float#

Return the cell volume.

Computed as the determinant of the basis-vector matrix. The unit follows that of the input cell (e.g., Angstrom^3 for an Angstrom-valued cell).

property Z: int#

Return the number of formula units in this cell.

Computed as the GCD of the per-species atom counts. For example, a cell with 4 Fe and 8 O atoms returns 4.

copy()[source]#

Return an independent copy of this cell.

Internal arrays (cell, scaled positions, masses, magnetic moments, species ids) are re-allocated by the constructor so that mutating one cell does not affect the other. The species table entries themselves are frozen dataclasses and are safe to share.

Return type:

PhonopyAtoms

totuple(distinguish_symbol_index=False)[source]#

Return (cell, scaled_positions, numbers).

If magnetic_moments is set, (cell, scaled_positions, numbers, magnetic_moments) is returned instead.

Parameters:

distinguish_symbol_index (bool, optional) – If True, the per-atom integer is the species id (atoms with the same symbol but different suffix get different ids); suitable as the types argument for spglib when symbol-suffix groupings must be preserved. If False (default), the per-atom integer is the atomic number. VCA cells always use species_ids regardless, since a VCA mixture has no single atomic number.

Returns:

(cell, scaled_positions, numbers) or (cell, scaled_positions, numbers, magnetic_moments).

Return type:

tuple

get_yaml_lines()[source]#

Return the crystal structure as a list of yaml text lines.

Return type:

list[str]

property formula: str#

Return chemical formula as a string.

The formula is constructed by sorting elements alphabetically and appending numbers for elements that appear more than once. E.g., “Si2O4” for two Si and four O atoms.

property reduced_formula: str#

Return reduced chemical formula as a string.

The reduced formula divides all element counts by their GCD. E.g., “Fe4O8” becomes “Fe2O4”.

property normalized_formula: str#

Return normalized formula as a string.

The normalized formula scales all element counts so they sum to 1. E.g., “Fe2O3” becomes “Fe0.4O0.6”.

Calculator interface helpers#

phonopy.interface.calculator.read_crystal_structure

Return crystal structure from file in each calculator format.

phonopy.interface.calculator.write_crystal_structure

Write crystal structure to file in each calculator format.

phonopy.interface.calculator.read_crystal_structure(filename=None, interface_mode=None, chemical_symbols=None, phonopy_yaml_cls=None)[source]#

Return crystal structure from file in each calculator format.

Parameters:
  • filename (str, optional) – Filename that contains cell structure information. Default is None. The predetermined filename for each interface_mode is used.

  • interface_mode (str, optional) – This is used to recognize the file format. Default is None, which is equivalent to ‘vasp’ mode.

  • chemical_symbols (list of str, optional) – This is only used for ‘vasp’ mode. VASP POSCAR file format can be written without chemical symbol information. With this option, chemical symbols can be given.

  • phonopy_yaml_cls (PhonopyYaml, optional) – This brings PhonopyYaml-like class dependent parameters. Here, currently only the default filenames are provided by this. Default is None.

Returns:

(Unit cell, structure info)

The StructureInfo contains unitcell_filename and, depending on the interface, additional calculator-specific fields.

Return type:

tuple[PhonopyAtoms | None, StructureInfo]

phonopy.interface.calculator.write_crystal_structure(filename, cell, interface_mode=None, optional_structure_info=None)[source]#

Write crystal structure to file in each calculator format.

filenamestr or os.PathLike

File name to be used to write out the crystal structure.

cellPhonopyAtoms

Crystal structure

interface_modestr, optional

Calculator interface such as ‘vasp’, ‘qe’, … Default is None, that is equivalent to ‘vasp’.

optional_structure_infoStructureInfo, optional

Information returned by the method read_crystal_structure. See the docstring. Default is None.

Parameters: