Equation of state#

The EquationOfState class can be used to find equilibrium volume, energy, and bulk modulus for solids:

class ase.eos.EquationOfState(volumes, energies, eos='sj')[source]#

Fit equation of state for bulk systems.

Parameters:

eos (str) –

Name of the equation of state to fit. Must be one of:

  • sj (default)

    A third order inverse polynomial fit 10.1103/PhysRevB.67.026103

                        2      3        -1/3
    E(V) = c + c t + c t  + c t ,  t = V
            0   1     2      3
    
  • taylor

    A third order Taylor series expansion about the minimum volume

  • murnaghan

    PRB 28, 5480 (1983)

  • birch

    Intermetallic compounds: Principles and Practice, Vol I: Principles. pages 195-210

  • birchmurnaghan

    PRB 70, 224107

  • pouriertarantola

    PRB 70, 224107

  • vinet

    PRB 70, 224107

  • anton-schmidt

    H. Anton and P. C. Schmidt, Intermetallics 5, 449 (1997). doi: 10.1016/s0966-9795(97)00017-4

  • p3

    A third order polynomial fit

Examples

>>> from ase.build import bulk
>>> from ase.calculators.emt import EMT
>>> atoms = bulk('Al')
>>> atoms.calc = EMT()
>>> cell = atoms.cell.copy()
>>> volumes = []
>>> energies = []
>>> for x in np.linspace(0.97, 1.01, 5):
...     atoms.set_cell(cell * x, scale_atoms=True)
...     volumes.append(atoms.get_volume())
...     energies.append(atoms.get_potential_energy())
>>> eos = EquationOfState(volumes, energies, eos='murnaghan')
>>> v0, e0, b0 = eos.fit()
>>> ax = eos.plot(show=False)
fit(warn=True)[source]#

Calculate volume, energy, and bulk modulus.

Returns the optimal volume, the minimum energy, and the bulk modulus. Notice that the ASE units for the bulk modulus is eV/Angstrom^3 - to get the value in GPa, do this:

v0, e0, B = eos.fit()
print(B / kJ * 1.0e24, 'GPa')
plot(filename=None, show=False, ax=None)[source]#

Plot fitted energy curve.

Uses Matplotlib to plot the energy curve. Use show=True to show the figure and filename=’abc.png’ or filename=’abc.eps’ to save the figure to a file.

Convenient helper function:

ase.eos.calculate_eos(atoms, npoints=5, eps=0.04, trajectory=None, callback=None)[source]#

Calculate equation-of-state.

atoms: Atoms object

System to calculate EOS for. Must have a calculator attached.

npoints: int

Number of points.

eps: float

Variation in volume from v0*(1-eps) to v0*(1+eps).

trajectory: Trjectory object or str

Write configurations to a trajectory file.

callback: function

Called after every energy calculation.

>>> from ase.build import bulk
>>> from ase.calculators.emt import EMT
>>> from ase.eos import calculate_eos
>>> a = bulk('Cu', 'fcc', a=3.6)
>>> a.calc = EMT()
>>> eos = calculate_eos(a, trajectory='Cu.traj')
>>> v, e, B = eos.fit()
>>> a = (4 * v)**(1 / 3.0)
>>> print('{0:.6f}'.format(a))
3.589825