Vibrational modes

You can calculate the vibrational modes of an Atoms object in the harmonic approximation using the Vibrations.

class ase.vibrations.Vibrations(atoms, indices=None, name='vib', delta=0.01, nfree=2)[source]

Class for calculating vibrational modes using finite difference.

The vibrational modes are calculated from a finite difference approximation of the Hessian matrix.

The summary(), get_energies() and get_frequencies() methods all take an optional method keyword. Use method=’Frederiksen’ to use the method described in:

T. Frederiksen, M. Paulsson, M. Brandbyge, A. P. Jauho: “Inelastic transport theory from first-principles: methodology and applications for nanoscale devices”, Phys. Rev. B 75, 205413 (2007)
atoms: Atoms object
The atoms to work on.
indices: list of int
List of indices of atoms to vibrate. Default behavior is to vibrate all atoms.
name: str
Name to use for files.
delta: float
Magnitude of displacements.
nfree: int
Number of displacements per atom and cartesian coordinate, 2 and 4 are supported. Default is 2 which will displace each atom +delta and -delta for each cartesian coordinate.

Example:

>>> from ase import Atoms
>>> from ase.calculators.emt import EMT
>>> from ase.optimize import BFGS
>>> from ase.vibrations import Vibrations
>>> n2 = Atoms('N2', [(0, 0, 0), (0, 0, 1.1)],
...            calculator=EMT())
>>> BFGS(n2).run(fmax=0.01)
BFGS:   0  16:01:21        0.440339       3.2518
BFGS:   1  16:01:21        0.271928       0.8211
BFGS:   2  16:01:21        0.263278       0.1994
BFGS:   3  16:01:21        0.262777       0.0088
>>> vib = Vibrations(n2)
>>> vib.run()
Writing vib.eq.pckl
Writing vib.0x-.pckl
Writing vib.0x+.pckl
Writing vib.0y-.pckl
Writing vib.0y+.pckl
Writing vib.0z-.pckl
Writing vib.0z+.pckl
Writing vib.1x-.pckl
Writing vib.1x+.pckl
Writing vib.1y-.pckl
Writing vib.1y+.pckl
Writing vib.1z-.pckl
Writing vib.1z+.pckl
>>> vib.summary()
---------------------
#    meV     cm^-1
---------------------
0    0.0       0.0
1    0.0       0.0
2    0.0       0.0
3    2.5      20.4
4    2.5      20.4
5  152.6    1230.8
---------------------
Zero-point energy: 0.079 eV
>>> vib.write_mode(-1)  # write last mode to trajectory file
clean(empty_files=False)[source]

Remove pickle-files.

Use empty_files=True to remove only empty files.

fold(frequencies, intensities, start=800.0, end=4000.0, npts=None, width=4.0, type='Gaussian', normalize=False)[source]

Fold frequencies and intensities within the given range and folding method (Gaussian/Lorentzian). The energy unit is cm^-1. normalize=True ensures the integral over the peaks to give the intensity.

get_energies(method='standard', direction='central', **kw)[source]

Get vibration energies in eV.

get_frequencies(method='standard', direction='central')[source]

Get vibration frequencies in cm^-1.

get_mode(n)[source]

Get mode number .

iterdisplace(inplace=False)[source]

Yield name and atoms object for initial and displaced structures.

Use this to export the structures for each single-point calculation to an external program instead of using run(). Then save the calculated gradients to <name>.pckl and continue using this instance.

iterimages()[source]

Yield initial and displaced structures.

run()[source]

Run the vibration calculations.

This will calculate the forces for 6 displacements per atom +/-x, +/-y, +/-z. Only those calculations that are not already done will be started. Be aware that an interrupted calculation may produce an empty file (ending with .pckl), which must be deleted before restarting the job. Otherwise the forces will not be calculated for that displacement.

Note that the calculations for the different displacements can be done simultaneously by several independent processes. This feature relies on the existence of files and the subsequent creation of the file in case it is not found.

If the program you want to use does not have a calculator in ASE, use iterdisplace to get all displaced structures and calculate the forces on your own.

summary(method='standard', direction='central', freq=None, log=<_io.StringIO object>)[source]

Print a summary of the vibrational frequencies.

Parameters:

method : string
Can be ‘standard’(default) or ‘Frederiksen’.
direction: string
Direction for finite differences. Can be one of ‘central’ (default), ‘forward’, ‘backward’.
freq : numpy array
Optional. Can be used to create a summary on a set of known frequencies.
log : if specified, write output to a different location than
stdout. Can be an object with a write() method or the name of a file to create.
write_dos(out='vib-dos.dat', start=800, end=4000, npts=None, width=10, type='Gaussian', method='standard', direction='central')[source]

Write out the vibrational density of states to file.

First column is the wavenumber in cm^-1, the second column the folded vibrational density of states. Start and end points, and width of the Gaussian/Lorentzian should be given in cm^-1.

write_jmol()[source]

Writes file for viewing of the modes with jmol.

write_mode(n=None, kT=0.02585199101165164, nimages=30)[source]

Write mode number n to trajectory file. If n is not specified, writes all non-zero modes.

name is a string that is prefixed to the names of all the files created. atoms is an Atoms object that is either at a fully relaxed ground state or at a saddle point. freeatoms is a list of atom indices for which the vibrational modes will be calculated, the rest of the atoms are considered frozen. displacements is a list of displacements, one for each free atom that are used in the finite difference method to calculate the Hessian matrix. method is -1 for backward differences, 0 for centered differences, and 1 for forward differences.