Modulation
Modulation#
- class spgrep_modulation.modulation.Modulation(primitive_symmetry: Symmetry, supercell: Supercell, dynamical_matrix: DynamicalMatrix | DynamicalMatrixNAC, qpoint: NDArrayFloat, nac_q_direction: NDArrayFloat | None = None, factor: float = 15.633302300230191, degeneracy_tolerance: float = 0.0001, seed: int = 0)[source]#
Generate modulated cells based on dynamical matrix of phonon.
Generate modulated cells based on dynamical matrix of phonon.
- Parameters
primitive_symmetry (phonopy.structure.symmetry.Symmetry) – phonopy’s Symmetry object for primitive cell
supercell (phonopy.structure.cells.Supercell) – phonopy’s supercell object to apply modulation
dynamical_matrix (phonopy.harmonic.dynamical_matrix.DynamicalMatrix) – phonopy’s dynamical matrix object
qpoint (array, (3, )) –
nac_q_direction – Used to calculate dynamical matrix around Gamma point if specified
factor (float) – Unit conversion
degeneracy_tolerance (float) – Absolute tolerance to groupby phonon frequencies in
factor
unitseed (int) – Seed value for sampling modulation from degenerated modes
- property dynamical_matrix: DynamicalMatrix | DynamicalMatrixNAC#
Return Phonopy’s dynamical matrix.
- property primitive_symmetry: phonopy.structure.symmetry.Symmetry#
Return Phonopy’s symmetry for primitive cell.
- property little_rotations: numpy.ndarray[Any, numpy.dtype[numpy.int32]]#
Return rotations of little group at
qpoint
.
- property little_translations: numpy.ndarray[Any, numpy.dtype[numpy.float64]]#
Return translations of little group at
qpoint
.
- property primitive: phonopy.structure.cells.Primitive#
Return Phonopy’s primitive cell.
- property supercell: phonopy.structure.cells.Supercell#
Return Phonopy’s supercell without modulations.
- property qpoint: numpy.ndarray[Any, numpy.dtype[numpy.float64]]#
Return qpoint.
- property unit_conversion_factor: float#
Return unit conversion for frequencies.
- property degeneracy_tolerance: float#
Return absolute tolerance to detect degenerated bands.
- property eigenspaces: list[tuple[float, numpy.ndarray[Any, numpy.dtype[numpy.complex128]], numpy.ndarray[Any, numpy.dtype[numpy.complex128]]]]#
Return list of
(eigenvalue, degenerated eigenvectors, irrep formed by the eigenvectors)
of dynamical matrix.
eigenvalue
: float
eigenvectors
: array, (dim, num_atoms, 3)
irrep
: array, (little_order, dim, dim)
- get_modulated_supercell_and_modulation(frequency_index: int, amplitudes: list[float], arguments: list[float], return_cell: bool = True) tuple[PhonopyAtoms, NDArrayComplex] | NDArrayComplex [source]#
Return modulated cell and modulation with given amplitudes and arguments.
- Parameters
frequency_index (int) – Index of frequency with dim eivenvectors
amplitudes (list[float], (dim, )) –
arguments (list[float], (dim, )) – in radian
return_cell (bool = True) –
- get_modulated_supercells(frequency_index: int, maximal_displacement: float = 0.11, max_size: int = 1) list[phonopy.structure.atoms.PhonopyAtoms] [source]#
Return modulated cells with randomly sampled arguments.
- Parameters
frequency_index (int) – Index of considered eigenmodes in
Modulation.eigenspaces
maximal_displacement (int) – Amplitude of modulation is chosen so that the maximal displacement in returned modulation is equal to
maximal_displacement
. Same unit asModulation.primitive.cell
.max_size (int) – Number of returned modulated structures when degenerated eigenmodes are specified.
- Returns
cells
- Return type
list of PhonopyAtoms
- get_high_symmetry_modulated_supercells(frequency_index: int, maximal_displacement: float = 0.11) list[phonopy.structure.atoms.PhonopyAtoms] [source]#
Return modulated cells corresponding to one-dimensional order-parameter directions of an isotropy subgroup.
- Parameters
frequency_index (int) – Index of considered eigenmodes in
Modulation.eigenspaces
maximal_displacement (int) – Amplitude of modulation is chosen so that the maximal displacement in returned modulation is equal to
maximal_displacement
. Same unit asModulation.primitive.cell
.- Returns
cells
- Return type
list of PhonopyAtoms
- classmethod with_supercell_and_symmetry_search(dynamical_matrix: DynamicalMatrix | DynamicalMatrixNAC, supercell_matrix: NDArrayInt, qpoint: NDArrayFloat, nac_q_direction: NDArrayFloat | None = None, factor: float = 15.633302300230191, degeneracy_tolerance: float = 0.0001, symprec: float = 1e-05, seed: int = 0) Modulation [source]#
Return Modulation class with supercell matrix.