Mode projection#

The mode projection functionality in dynasor is mainly handled by two objects: The ModeProjector class and project_modes() function. The ModeProjector provides access to data objects representing a q-point QPoint and from the q-point there is access to an object representing a particular band Band at that q-point. In addition, simple wrappers around the coordinates Q, P and F exist via ComplexCoordinate to easily set the amplitude and phase of a mode while preserving the \(Q(-q)=Q^*(q)\) symmetries. Internally dynasor wraps the primitive cell Prim and supercell Supercell. As a user only the ModeProjector and project_modes should need to be imported.

class dynasor.ModeProjector(primitive, supercell, force_constants)[source]#

The ModeProjector maps between real atomic displacements u and complex mode coordinates Q.

Some special python methods are implemented. The __str__ and __repr__ provides useful info. QPoint objects are representations of a single q-point and associated information and can be accessed either by call providing a reduced wavevector

>>> mp((1/2, 0, 0))  

or by index corresponding to reduced q-point accessible from q_reduced

>>> mp[2]  

In addition to mode coordinates Q the class can also map the atomic velocities v to mode momenta P as well as atomic forces f to mode forces F. The class can also map back. This mapping is done using getters and setters. Internally only Q, P and F are stored

>>> Q = mp.get_Q()  
>>> mp.set_u(u)  

In addition, the forces corresponding to the harmonic forces can be accessed by get_f_harmonic() and get_F_harmonic(). For ASE Atoms objects the displacments etc. can be updated and applied by

>>> mp.update_from_atoms(atoms)  
>>> atoms = mp.get_atoms(harmonic_forces=False)  

The shapes for each property uses the follwoing varaibles

  • N: Number of atoms in supercell

  • Nu: unit cells (N/Np)

  • Np: primitive basis atoms (N/Nu)

  • Nb: bands (Np*3)

  • Nq: q-points (Nu)

Please consult the documentation or the specific getters and setters to see the exact transformations used.

Units#

The internal units in dynasor are Å, fs and eV. All frequencies are angular (a.k.a. the “physicist’s convention” with 2π included). These are the units dynasor will expect and return. In, e.g., print functions conventional units such fs, Å, THz, Da, meV are commonly used.

Mass: The internal unit choice (eV, Å, fs) means that the mass unit is not Dalton but rather 0.009648533290731906 Da. We refer to this unit as the “dynasor mass unit” (dmu), i.e., 1 Da = 103.64269572045423 dmu. As a user you will only see this unit in the output of ModeProjector objects. Masses provided via, e.g., ASE Atoms objects are converted internally.

Waves: dynasor reports and expects spatial (angular) frequencies in rad/Å and temporal (angular) frequencies in rad/fs. This follows the often-used convention in physics to include the factor of 2π in the wave vectors. For instance the wavelength is given by λ=2π/q.

Mode amplitudes: Mode amplitudes are reported in Å√dmu = fs√eV.

Velocities: For modes the momenta are reported in Å√dmu/fs or just √eV while atomic velocities are reported in Å/fs.

Mode forces: The force is defined as the derivative of the momenta with respect to time so the unit used when reporting mode forces is Å√dmu/fs² (or √eV/fs).

Internal arrays#

For the curious, the internal data arrays are

  • primitive, supercell, force_constants (input)

  • _q, q_minus (reduced q-points and which q-points are related by inversion)

  • _D, _w2, _W (dynamical matrices, frequencies (ev/Ų/Da), polarization vectors)

  • _X (eigenmodes which are mass weighted “polarization vectors” in the supercell)

property eigenmodes: ndarray[tuple[int, ...], dtype[float]]#

The eigenmodes in the supercell as (Nq, Nb, N, 3)-array

The eigenmodes include the masses such that \(Q = X u\) where \(u\) are the supercell displacements.

get_F()[source]#

The mode force amplitudes in eV/Å√dmu.

Return type:

ndarray[tuple[int, ...], dtype[float]]

get_F_harmonic()[source]#

The harmonic mode forces, computed as \(-\omega^2 * Q\) in eV/Å√dmu.

Return type:

ndarray[tuple[int, ...], dtype[float]]

get_P()[source]#

The mode momentum amplitudes in √eV.

Return type:

ndarray[tuple[int, ...], dtype[float]]

get_Q()[source]#

The mode coordinate amplitudes in Å√dmu.

Return type:

ndarray[tuple[int, ...], dtype[float]]

get_atoms(harmonic_forces=False)[source]#

Returns ASE Atoms object with displacement, velocities, forces, and harmonic energies.

Parameters:

harmonic_forces (Optional[bool]) – Whether the forces should be taken from the internal F or via -omega^2 Q.

Return type:

Atoms

get_f()[source]#

The atomic forces in eV/Å.

Return type:

ndarray[tuple[int, ...], dtype[float]]

get_f_harmonic()[source]#

The harmonic atomic forces for the current displacements.

Return type:

ndarray[tuple[int, ...], dtype[float]]

get_u()[source]#

The atomic displacements in Å.

Return type:

ndarray[tuple[int, ...], dtype[float]]

get_v()[source]#

The atomic velocities in Å/fs.

Return type:

ndarray[tuple[int, ...], dtype[float]]

property kinetic_energies: ndarray[tuple[int, ...], dtype[float]]#

Kinetic energy per mode as (Nq, Nb)-array.

The kinetic energies are defined as \(1/2 P P^*\). Should equal \(1/2 k_B T\) in equilibrium.`

property omegas: ndarray[tuple[int, ...], dtype[float]]#

The frequencies of each mode in rad/fs.

Following convetion negative values indicate imaginary frequencies.

property polarizations: ndarray[tuple[int, ...], dtype[float]]#

The polarization vectors for each mode (Nq, Nb, Np, 3).

property potential_energies: ndarray[tuple[int, ...], dtype[float]]#

Potential energy per mode as (Nq, Nb)-array.

The potential energies are defined as \(1/2 Q Q^*\) and should equal \(1/2 k_B T\) in equilibrium for a harmonic system.

property q_cartesian: ndarray[tuple[int, ...], dtype[float]]#

The q-points in cartesian coordinates with unit of rad/Å (2π included).

property q_minus: ndarray[tuple[int, ...], dtype[float]]#

The index of the corresponding counter-propagating mode (\(-q\)).

property q_reduced: ndarray[tuple[int, ...], dtype[float]]#

The q-points in reduced coordinates.

For example a zone boundary mode would be (1/2, 0, 0)

classmethod read(file_name)[source]#

Return ModeProjector instance from pickle file that was saved using write().

Return type:

ModeProjector

set_F(F)[source]#

Sets the internal mode forces \(F\).

The functions ensures \(F(-q)=F^*(q)\).

Return type:

None

set_P(P)[source]#

Sets the internal mode momenta \(P\).

The functions ensures \(P(-q)=P^*(q)\)

Return type:

None

set_Q(Q)[source]#

Sets the internal mode coordinates \(Q\).

The functions ensures \(Q(-q)=Q^*(q)\).

Return type:

None

set_f(f)[source]#

Sets the internal mode forces \(F\) given the atomic forces \(f\).

\[F = X^* * f / m\]
Parameters:

f (ndarray[tuple[int, ...], dtype[float]]) – The atomic forces in eV/Å.

Return type:

None

set_u(u)[source]#

Sets the internal mode coordinates \(Q\) given the atomic displacements \(u\).

\[Q = X u\]
Parameters:

u (ndarray[tuple[int, ...], dtype[float]]) – The atomic displacements in Å.

Return type:

None

set_v(v)[source]#

Sets the internal mode momenta \(P\) given the atomic velocities \(v\).

\[P = X^* * v\]
Parameters:

v (ndarray[tuple[int, ...], dtype[float]]) – The atomic velocities in Å/fs.

Return type:

None

update_from_atoms(atoms)[source]#

Updates the ModeProjector objects with displacments, velocities, and forces from an ASE Atoms object.

Checks for an attached calculator in the first place and next for a forces array.

If no data sets corresponding array to zeros.

The masses and velocities are converted to dynasor units internally.

Return type:

None

property virial_energies: ndarray[tuple[int, ...], dtype[float]]#

The virial energies per mode as (Nq, Nb)-array.

The virial energies are defined here as \(-1/2 Q F\), which should have an expectation value of \(1/2 k_B T\) per mode in equilibrium. For a harmonic system this is simply equal to the potential energy. This means that the virial energy can be used to monitor the anharmonicity or define a measure of the potential energy.

write(file_name)[source]#

Uses pickle to write mode projector to file.

Return type:

None

class dynasor.modes.qpoint.QPoint(q_index, mp)[source]#

Representation of a single q-point and properties.

The bands can be accessed by, e.g., qp[2] to get band index 2 in the form of a Band object.

Parameters:
  • q – q-point index.

  • mp – Mode project object.

property eigenmodes: ndarray[tuple[int, ...], dtype[float]]#

Slice, see ModeProjector.

property index: int#

q-point index corresponding to q_reduced.

property is_real: bool#

If the q-point has purely real mode coordinates, q=-q.

property kinetic_energies: ndarray[tuple[int, ...], dtype[float]]#

Slice, see ModeProjector.

property omegas: ndarray[tuple[int, ...], dtype[float]]#

Slice, see ModeProjector.

property polarizations: ndarray[tuple[int, ...], dtype[float]]#

Slice, see ModeProjector.

property potential_energies: ndarray[tuple[int, ...], dtype[float]]#

Slice, see ModeProjector.

property q_cartesian: ndarray[tuple[int, ...], dtype[float]]#

Slice, see ModeProjector.

property q_minus#

The corresponding counter-propagating mode.

property q_reduced: ndarray[tuple[int, ...], dtype[float]]#

Slice, see ModeProjector.

property virial_energies: ndarray[tuple[int, ...], dtype[float]]#

Slice, see ModeProjector.

property wavelength: float#

Wavelength of mode in Å.

property wavenumber: float#

Wavenumber of mode in rad/Å.

class dynasor.modes.band.Band(q_index, band_index, mp)[source]#

Represents properties of a single band belonging to a q-point.

To modify the coordinates Q, P and F use band.Q and see doc for ComplexCoordinate.

Parameters:
  • q_index (int) – q-point index.

  • band_index (int) – Band index.

  • mp (ModeProjector) – Mode projector.

property F: ndarray[tuple[int, ...], dtype[float]]#

The mode force

property P: float#

The mode momentum.

property Q: float#

The mode coordinate.

property eigenmode: ndarray[tuple[int, ...], dtype[float]]#

Slice, see doc for ModeProjector.

property index: int#

The index of the band at the q-point.

property kinetic_energy: float#

Slice, see doc for ModeProjector.

property omega: float#

Slice, see doc for ModeProjector.

property polarization: ndarray[tuple[int, ...], dtype[float]]#

Slice, see doc for ModeProjector.

property potential_energy: float#

Slice, see doc for ModeProjector.

property qpoint: QPoint#

The q-point to which the band belongs.

property virial_energy: float#

Slice, see doc for ModeProjector

class dynasor.modes.complex_coordinate.ComplexCoordinate(q, s, mp)[source]#

Class to work with the in general complex coordinates of lattice dynamics.

Can be cast to complex by complex(cc)

Example

>>> cc.complex = 1.7j  
>>> cc.amplitude = 3.8  
>>> cc.phase = 2*np.pi  
property amplitude#

The amplitude of the complex coordinate.

property band: Band#

The band to which this coordinate belongs to.

property complex#

The complex coordinate as a complex number.

property phase#

The phase of the complex coordinate.

dynasor.project_modes(traj, modes, ideal_supercell, check_mic=True)[source]#

Projects an atomic trajectory onto set of phonon modes.

Parameters:
  • traj (Trajectory) – Input trajectory.

  • modes (ndarray[tuple[int, ...], dtype[float]]) – Modes to project on (1, Nlambda, N, 3) array where Nm is the number of modes and N is the number of atoms in the supercell.

  • ideal_supercell (Atoms) – Ideal supercell used to find atomic displacements and should correspond to the ideal structure. Be careful not to mess up the permutation.

  • check_mic (Optional[bool]) – Whether to wrap the displacements or not, faster if no wrap.

Return type:

tuple[ndarray[tuple[int, ...], dtype[float]], ndarray[tuple[int, ...], dtype[float]]]

Returns:

A tuple comprising (Q,P) where Q are the mode coordinates as a complex array with dimension (length of traj, number of modes) and P are the mode momenta as a complex array with dimension (length of traj, number of modes).

class dynasor.modes.atoms.Prim(atoms)[source]#
class dynasor.modes.atoms.Supercell(supercell, prim)[source]#

The supercell takes care of some mappings between the primitive and repeated structure.

In particular the P-matrix connecting the cells as well as the offset-index of each atom is calculated.

Note that the positions cannot be revovered as offset x cell + basis since the atoms get wrapped.

Parameters:
  • supercell (Union[Atoms, DynasorAtoms]) – Some ideal repetition of the primitive structure and possible wrapping.

  • prim (Union[Atoms, DynasorAtoms]) – Primitive structure.

property P: ndarray[tuple[int, ...], dtype[float]]#

P-matrix is defined as dot(P, prim.cell) = supercell.cell

property P_inv: ndarray[tuple[int, ...], dtype[float]]#

Inverse of P.

property indices: ndarray[tuple[int, ...], dtype[int]]#

The basis index of each atom

property n_cells: int#

Number of unit cells

property offsets: ndarray[tuple[int, ...], dtype[float]]#

The offset of each atom.