Modulation¶
- class spgrep_modulation.modulation.Modulation(primitive_symmetry, supercell, dynamical_matrix, qpoint, nac_q_direction=None, factor=15.633302300230191, degeneracy_tolerance=0.0001, seed=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¶
Return Phonopy’s dynamical matrix.
- property primitive_symmetry¶
Return Phonopy’s symmetry for primitive cell.
- property little_rotations¶
Return rotations of little group at
qpoint
.
- property little_translations¶
Return translations of little group at
qpoint
.
- property primitive¶
Return Phonopy’s primitive cell.
- property supercell¶
Return Phonopy’s supercell without modulations.
- property qpoint¶
Return qpoint.
- property unit_conversion_factor¶
Return unit conversion for frequencies.
- property degeneracy_tolerance¶
Return absolute tolerance to detect degenerated bands.
- property eigenspaces¶
Return list of
(eigenvalue, degenerated eigenvectors, irrep formed by the eigenvectors)
of dynamical matrix.
frequency_index
is the index ofeigenspaces
. Eachfrequency_index
corresponds to a eigenvalue of the dynamical matrix, \(\omega_{\alpha s}^{2}\). Here, \(\alpha\) specifies irrep up to unitary equivalence and \(s\) distinguishes the irrep among equivalent irreps.eigenspaces[frequency_index]
is a tuple with
eigenvalue
: float, eigenvalue of the dynamical matrix \(\omega_{\alpha s}^{2}\)
eigenvectors
: array, (dim, num_atoms, 3),eigenvectors[lambda, kappa, :]
is an eigenvector \(\mathbf{e}(\kappa; \mathbf{q}\alpha s \lambda)\).
irrep
: array, (little_order, dim, dim), representation matrices of the irrep \(\tilde{\Gamma}^{\mathbf{q}\alpha s}\)
- get_modulated_supercell_and_modulation(frequency_index, amplitudes, arguments, return_cell=True)[source]¶
Return modulated cell and modulation with given amplitudes and arguments.
- Parameters
frequency_index (int) – Index of considered eigenmodes in
Modulation.eigenspaces
. SeeModulation.eigenspaces()
.amplitudes (list[float], (dim, )) –
arguments (list[float], (dim, )) – in radian
return_cell (bool = True) –
- get_modulated_supercells(frequency_index, maximal_displacement=0.11, max_size=1)[source]¶
Return modulated cells with randomly sampled arguments.
- Parameters
frequency_index (int) – Index of considered eigenmodes in
Modulation.eigenspaces
. SeeModulation.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, maximal_displacement=0.11)[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
. SeeModulation.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