Analysis tools

Thermodynamics

Formation Energies

class qmpy.FormationEnergy(*args, **kwargs)[source]

Base class for a formation energy.

Relationships:
Calculation via calculation
Composition via composition
Entry via entry
FormationEnergy via equilibrium
Fit via fit
Attributes:
id
delta_e: Formation energy (eV/atom)
description: A label of the source of the formation energy.
stability: Distance from the convex hull (eV/atom)
exception DoesNotExist
exception MultipleObjectsReturned
save(*args, **kwargs)[source]

Save the current instance. Override this in a subclass if you want to control the saving process.

The ‘force_insert’ and ‘force_update’ parameters can be used to insist that the “save” must be an SQL insert or update (or equivalent for non-SQL backends), respectively. Normally, they should not be set.

class qmpy.ExptFormationEnergy(*args, **kwargs)[source]

Experimentally measured formation energy.

Any external formation energy should be entered as an ExptFormationEnergy object, rather than a FormationEnergy. If the external source is also computational, set the “dft” attribute to be True.

Relationships:
Composition via composition
Fit via fit
Attributes:
id: integer primary key.
delta_e: measured formation energy.
delta_g: measured free energy of formation.
dft: (bool) True if the formation energy is from a non-OQMD DFT
calculation.
source: (str) Identifier for the source.
exception DoesNotExist
exception MultipleObjectsReturned

Reference energies

class qmpy.Fit(*args, **kwargs)[source]

The core model for a reference energy fitting scheme.

The Fit model links to the experimental data (ExptFormationEnergy objects) that informed the fit, as well as the DFT calculations (Calculation objects) that were matched to each experimental formation energy. Once the fit is completed, it also stores a list of chemical potentials both as a relationship to ReferenceEnergy and HubbardCorrection objects. These correction energies can also be accessed by dictionaries at Fit.mus and Fit.hubbard_mus.

Relationships:
Calculation via dft
ExptFormationEnergy via experiments
FormationEnergy via formationenergy_set
HubbardCorrection via hubbard_correction_set
ReferenceEnergy via reference_energy_set
Attributes:
name: Name for the fitting

Examples:

>>> f = Fit.get('standard')
>>> f.experiments.count()
>>> f.dft.count()
>>> f.mus
>>> f.hubbard_mus
exception DoesNotExist
exception MultipleObjectsReturned
class qmpy.ReferenceEnergy(*args, **kwargs)[source]

Elemental reference energy for evaluating heats of formation.

Relationships:
Fit via fit
Element via element
Attributes:
id
value: Reference energy (eV/atom)
exception DoesNotExist
exception MultipleObjectsReturned
class qmpy.HubbardCorrection(*args, **kwargs)[source]

Energy correction for DFT+U energies.

Relationships:
Fit via fit
Element via element
Hubbard via hubbard
Attributes:
id
value: Correction energy (eV/atom)
exception DoesNotExist
exception MultipleObjectsReturned

Phase Space

class qmpy.PhaseSpace(bounds, mus=None, data=None, **kwargs)[source]

A PhaseSpace object represents, naturally, a region of phase space.

The most fundamental property of a PhaseSpace is its bounds, which are given as a hyphen-delimited list of compositions. These represent the extent of the phase space, and determine which phases are within the space.

Next, a PhaseSpace has an attribute, data, which is a PhaseData object, and is a container for Phase objects, which are used when performing thermodynamic analysis on this space.

The majority of attributes are lazy, that is, they are only computed when they are requested, and how to get them (of which there are often several ways) is decided based on the size and shape of the phase space.

property bound_elements

Alphabetically ordered list of elements with constrained composition.

property bound_space

Set of elements _of fixed composition_ in the PhaseSpace.

Examples:

>>> s = PhaseSpace('Fe-Li', 'O=-1.4')
>>> s.bound_space
set(['Fe', 'Li'])
property chempot_dimension

Chemical potential dimension.

Examples:

>>> s = PhaseSpace('Fe-Li', 'O=-2.5')
>>> s.chempot_dimension
0
>>> s = PhaseSpace('Fe-Li', 'N=0:-5')
>>> s.chempot_dimension
1
>>> s = PhaseSpace('Fe-Li', 'N=0:-5 F=0:-5')
>>> s.chempot_dimension
2
chempot_scan(element=None, umin=None, umax=None)[source]

Scan through chemical potentials of element from umin to umax identifing values at which phase transformations occur.

clear_all()[source]

Clears input data and analyzed results. Same as: >>> PhaseData.clear_data() >>> PhaseData.clear_analysis()

clear_analysis()[source]

Clears all calculated results.

clear_data()[source]

Clears all phase data.

property cliques

Iterator over maximal cliques in the phase space. To get a list of cliques, use list(PhaseSpace.cliques).

comp(coord)[source]

Returns the composition of a coordinate in phase space.

Examples:

>>> space = PhaseSpace('Fe-Li-O')
>>> space.comp([0.2, 0.2, 0.6])
{'Fe': 0.2, 'O': 0.6, 'Li': 0.2}
property comp_dimension

Compositional dimension of the region of phase space.

Examples:

>>> s = PhaseSpace('Fe-Li-O')
>>> s.comp_dimension
2
>>> s = PhaseSpace('FeO-Ni2O-CoO-Ti3O4')
>>> s.comp_dimension
3
compute_formation_energies()[source]

Evaluates the formation energy of every phase with respect to the chemical potentials in the PhaseSpace.

compute_stabilities(phases=None, save=False, reevaluate=True)[source]

Calculate the stability for every Phase.

Keyword Arguments:
phases:

List of Phases. If None, uses every Phase in PhaseSpace.phases

save:

If True, save the value for stability to the database.

new_only:

If True, only compute the stability for Phases which did not import a stability from the OQMD. False by default.

compute_stability(p)[source]

Compute the energy difference between the formation energy of a Phase, and the energy of the convex hull in the absence of that phase.

coord(composition, tol=0.0001)[source]

Returns the barycentric coordinate of a composition, relative to the bounds of the PhaseSpace. If the object isn’t within the bounds, raises a PhaseSpaceError.

Examples:

>>> space = PhaseSpace('Fe-Li-O')
>>> space.coord({'Fe':1, 'Li':1, 'O':2})
array([ 0.25,  0.25,  0.5 ])
>>> space = PhaseSpace('Fe2O3-Li2O')
>>> space.coord('Li5FeO4')
array([ 0.25,  0.75])
property dual_spaces

List of sets of elements, such that any possible tie-line between two phases in phases is contained in at least one set, and no set is a subset of any other.

property elements

Alphabetically ordered list of elements present in the PhaseSpace.

find_reaction_mus(element=None)[source]

Find the chemical potentials of a specified element at which reactions occur.

Examples:

>>> s = PhaseSpace('Fe-Li-O')
>>> s.find_reaction_mus('O')
gclp(composition={}, mus={}, phases=[])[source]

Returns energy, phase composition which is stable at given composition

Examples:

>>> space = PhaseSpace('Fe-Li-O')
>>> energy, phases = space.gclp('FeLiO2')
>>> print phases
>>> print energy
get_hull_points()[source]

Gets out-of PhaseSpace points. i.e. for FeSi2-Li, there are no other phases in the space, but there are combinations of Li-Si phases and Fe-Si phases. This method returns a list of phases including composite phases from out of the space.

Examples:

>>> space = PhaseSpace('FeSi2-Li')
>>> space.get_hull_points()
[<Phase FeSi2 (23408): -0.45110217625>,
<Phase Li (104737): 0>,
<Phase 0.680 Li13Si4 + 0.320 FeSi : -0.3370691816>,
<Phase 0.647 Li8Si3 + 0.353 FeSi : -0.355992801765>,
<Phase 0.133 Fe3Si + 0.867 Li21Si5 : -0.239436904167>,
<Phase 0.278 FeSi + 0.722 Li21Si5 : -0.306877209723>]
get_minima(phases, bounds)[source]

Given a set of Phases, get_minima will determine the minimum free energy elemental composition as a weighted sum of these compounds

get_phase_diagram(**kwargs)[source]

Creates a Renderer attribute with appropriate phase diagram components.

Examples:

>>> space = PhaseSpace('Fe-Li-O')
>>> space.get_renderer()
>>> plt.show()
get_qhull(phases=None, mus={})[source]

Get the convex hull for a given space.

get_reaction(var, facet=None)[source]

For a given composition, what is the maximum delta_composition reaction on the given facet. If None, returns the whole reaction for the given PhaseSpace.

Examples:

>>> space = PhaseSpace('Fe2O3-Li2O')
>>> equilibria = space.hull[0]
>>> space.get_reaction('Li2O', facet=equilibria)
get_reactions(var, electrons=1.0)[source]

Returns a list of Reactions.

Examples:

>>> space = PhaseSpace('Fe-Li-O')
>>> space.get_reactions('Li', electrons=1)
get_tie_lines_by_gclp(iterable=False)[source]

Runs over pairs of Phases and tests for equilibrium by GCLP. Not recommended, it is very slow.

property graph

networkx.Graph representation of the phase space.

property hull

List of facets of the convex hull.

in_bounds(composition)[source]

Returns True, if the composition is within the bounds of the phase space

Examples:

>>> space = PhaseSpace('Fe2O3-NiO2-Li2O')
>>> space.in_bounds('Fe3O4')
False
>>> space.in_bounds('Li5FeO8')
True
in_space(composition)[source]

Returns True, if the composition is in the right elemental-space for this PhaseSpace.

Examples:

>>> space = PhaseSpace('Fe-Li-O')
>>> space.in_space('LiNiO2')
False
>>> space.in_space('Fe2O3')
True
load(**kwargs)[source]

Loads oqmd data into the associated PhaseData object.

make_1d_vs_chempot(**kwargs)[source]

Plot of phase stability vs chemical potential for a single composition.

Examples:

>>> s = PhaseSpace('Fe', mus={'O':[0,-4]})
>>> r = s.make_vs_chempot()
>>> r.plot_in_matplotlib()
>>> plt.show()
make_as_binary(**kwargs)[source]

Construct a binary phase diagram (convex hull) and write it to a Renderer.

Examples:

>>> s = PhaseSpace('Fe-P')
>>> r = s.make_as_binary()
>>> r.plot_in_matplotlib()
>>> plt.show()
make_as_graph(**kwargs)[source]

Construct a graph-style visualization of the phase diagram.

make_as_quaternary(**kwargs)[source]

Construct a quaternary phase diagram and write it to a Renderer.

Examples:

>>> s = PhaseSpace('Fe-Li-O-P')
>>> r = s.make_as_quaternary()
>>> r.plot_in_matplotlib()
>>> plt.show()
make_as_ternary(**kwargs)[source]

Construct a ternary phase diagram and write it to a Renderer.

Examples:

>>> s = PhaseSpace('Fe-Li-O-P')
>>> r = s.make_as_quaternary()
>>> r.plot_in_matplotlib()
>>> plt.show()
make_as_unary(**kwargs)[source]

Plot of phase volume vs formation energy.

Examples:

>>> s = PhaseSpace('Fe2O3')
>>> r = s.make_as_unary()
>>> r.plot_in_matplotlib()
>>> plt.show()
make_vs_chempot(**kwargs)[source]

Plot of phase stability vs chemical potential for a range of compositions.

Examples:

>>> s = PhaseSpace('Fe-Li', mus={'O':[0,-4]})
>>> r = s.make_vs_chempot()
>>> r.plot_in_matplotlib()
>>> plt.show()
property phase_diagram

Renderer of a phase diagram of the PhaseSpace

plot_reactions(var, electrons=1.0, save=False)[source]

Plot the convex hull along the reaction path, as well as the voltage profile.

save_tie_lines()[source]

Save all tie lines in this PhaseSpace to the OQMD. Stored in Formation.equilibrium

property shape

(# of compositional dimensions, # of chemical potential dimensions) The shape attribute of the PhaseSpace determines what type of phase diagram will be drawn.

Examples:

>>> s = PhaseSpace('Fe-Li', 'O=-1.2')
>>> s.shape
(1, 0)
>>> s = PhaseSpace('Fe-Li', 'O=0:-5')
>>> s.shape
(1, 1)
>>> s = PhaseSpace('Fe-Li-P', 'O=0:-5')
>>> s.shape
(2,1)
>>> s = PhaseSpace('Fe', 'O=0:-5')
>>> s.shape
(0, 1)
property space

Set of elements present in the PhaseSpace.

Examples:

>>> s = PhaseSpace('Pb-Te-Se')
>>> s.space
set(['Pb', 'Te', 'Se'])
>>> s = PhaseSpace('PbTe-Na-PbSe')
>>> s.space
set(['Pb', 'Te', 'Na', 'Se'])
property spaces

List of lists of elements, such that every phase in self.phases is contained in at least one set, and no set is a subset of any other. This corresponds to the smallest subset of spaces that must be analyzed to determine the stability of every phase in your dataset.

Examples:

>>> pa, pb, pc = Phase('A', 0), Phase('B', 0), Phase('C', 0)
>>> p1 = Phase('AB2', -1)
>>> p2 = Phase('B3C', -4)
>>> s = PhaseSpace('A-B-C', load=None)
>>> s.phases = [ pa, pb, pc, p1, p2 ]
>>> s.spaces
[['C', 'B'], ['A', 'B']]
stability_range(p, element=None)[source]

Calculate the range of phase p with respect to element.

property stable

List of stable phases

property tie_lines

List of length 2 tuples of phases with tie lines between them

property unstable

List of unstable phases.

class qmpy.PhaseData[source]

A PhaseData object is a container for storing and organizing phase data. Most importantly used when doing a large number of thermodynamic analyses and it is undesirable to access the database for every space you want to consider.

add_phase(phase)[source]

Add a phase to the PhaseData collection. Updates the PhaseData.phase_dict and PhaseData.phases_by_elt dictionaries appropriately to enable quick access.

Examples:

>>> pd = PhaseData()
>>> pd.add_phase(Phase(composition='Fe2O3', energy=-3))
>>> pd.add_phase(Phase(composition='Fe2O3', energy=-4))
>>> pd.add_phase(Phase(composition='Fe2O3', energy=-5))
>>> pd.phase_dict
{'Fe2O3': <Phase Fe2O3 : -5}
>>> pd.phases_by_elt['Fe']
[<Phase Fe2O3 : -3>, <Phase Fe2O3 : -4>, <Phase Fe2O3 : -5>]
add_phases(phases)[source]

Loops over a sequence of phases, and applies add_phase to each.

Equivalent to:

>>> pd = PhaseData()
>>> for p in phases:
>>>     pd.add_phase(p)
dump(filename=None, minimal=True)[source]

Writes the contents of the phase data to a file or to stdout.

Keyword Arguments:
filename:

If None, prints the file to stdout, otherwise writes the file to the specified filename. Default=None.

minimal:

Dump _every_ phase in the PhaseData object, or only those that can contribute to a phase diagram. If True, only the lowest energy phase at a given composition will be written. Default=True.

get_phase_data(space)[source]

Using an existing PhaseData object return a PhaseData object which is populated by returning a subset which is inside a given region of phase space.

Arguments:

space: formatted as in qmpy.PhaseSpace.__init__()

Examples:

>>> pd = PhaseData()
>>> pd.read_file('legacy.dat')
>>> new_pd = pd.get_phase_data(['Fe', 'O'])
>>> new_pd.phase_dict
load_library(library)[source]

Load a library file, containing self-consistent thermochemical data.

Equivalent to::
>>> pd = PhaseData()
>>> pd.read_file(INSTALL_PATH+'/data/thermodata/%s' % library)
load_oqmd(space=None, search={}, exclude={}, stable=False, fit='standard', total=False)[source]

Load data from the OQMD.

Keyword Arguments:
space:

sequence of elements. If supplied, will return only phases within that region of phase space. i.e. [‘Fe’, ‘O’] will return Fe, O and all iron oxides.

search:

dictionary of database search keyword:value pairs.

stable:

Restrict search to only stable phases (faster, but relies on having current phase stability analyses).

Examples::
>>> pd = PhaseData()
>>> search = {'calculation__path__contains':'icsd'}
>>> pd.load_oqmd(space=['Fe','O'], search=search, stable=True)
property phases

List of all phases.

read_file(filename, per_atom=True)[source]

Read in a thermodata file (named filename).

File format:

composition energy
Fe 0.0
O 0.0
Li 0.0
Fe3O4 -0.21331204979
FeO -0.589343204057
Fe3O4 -0.21331204979
FeLiO2 -0.446739168889
FeLi5O4 -0.198830531099
Keyword Arguments:

per_atom: If True, the supplied energies are per atom, not per formula unit. Defaults to True.

Symmetry

Spacegroup

class qmpy.Spacegroup(*args, **kwargs)[source]

Base class for a space group.

Relationships:
Structure via structure_set
Translation via centering_vectors
Operation via operations
WyckoffSite via site_set
Attributes:
number: Spacegroup #. (primary key)
centrosymmetric: (bool) Is the spacegroup centrosymmetric.
hall: Hall symbol.
hm: Hermann-Mauguin symobl.
lattice_system: Cubic, Hexagonal, Tetragonal, Orthorhombic,
Monoclinic or Triclinic.
pearson: Pearson symbol
schoenflies: Schoenflies symbol.
exception DoesNotExist
exception MultipleObjectsReturned
get_site(symbol)[source]

Gets WyckoffSite by symbol.

property rotations

List of rotation operations for the spacegroup.

save(*args, **kwargs)[source]

Save the current instance. Override this in a subclass if you want to control the saving process.

The ‘force_insert’ and ‘force_update’ parameters can be used to insist that the “save” must be an SQL insert or update (or equivalent for non-SQL backends), respectively. Normally, they should not be set.

property sym_ops

List of (rotation, translation) pairs for the spacegroup

property symbol

Returns the Hermann-Mauguin symbol for the spacegroup

property translations

List of translation operations for the spacegroup.

property wyckoff_sites

List of WyckoffSites.

Wyckoff Site

class qmpy.WyckoffSite(*args, **kwargs)[source]

Base class for a Wyckoff site. (e.g. a “b” site).

Relationships:
Spacegroup via spacegroup
Atom via atom_set
Site via site_set
Attributes:
id
symbol: Site symbol
multiplicity: Site multiplicity
x, y, z: Coordinate symbols.
exception DoesNotExist
exception MultipleObjectsReturned

Symmetry Operations

class qmpy.Operation(*args, **kwargs)[source]

A symmetry operation (rotation + translation).

Relationships:
Spacegroup via spacegroup
Rotation via rotation_set
Translation via translation_set
Attributes:
id

Examples:

>>> op = Operation.get('x+y-1/2,-z-y+1/2,x-z+1/2')
>>> print op
<Operation: +x+y+1/2,-y-z+1/2,+x-z+1/2>
exception DoesNotExist
exception MultipleObjectsReturned
classmethod get(value)[source]

Accepts symmetry operation strings, i.e. “+x, x+1/2, x+y-z” or a tuple of rotation matrix and translation vector.

Example:

>>> Operation.get("x,y,-y")
>>> Operation.get(( rot, trans ))
class qmpy.Translation(*args, **kwargs)[source]

A translation operation.

Relationships:
Spacegroup via spacegroup
Operation via operation
Attributes:
id
x, y, z: Translation vector. Accessed via vector.

Examples:

>>> trans = Translation.get([0, 0, 0.5])
<Translation: 0,0,+1/2>
>>> print(trans)
0,0,+1/2
>>> print(trans.vector)
array([ 0. ,  0. ,  0.5])
exception DoesNotExist
exception MultipleObjectsReturned
class qmpy.Rotation(*args, **kwargs)[source]

A rotation operation.

Relationships:
Spacegroup via spacegroup
Operation via operation
Attributes:
id
a11, a12, a13
a21, a22, a23
a31, a32, a33: Rotation matrix. Accessed via matrix.

Examples:

>>> rot = Rotation.get([[1, 0, 0], [-1, -1, 0], [0, 1, -1]])
<Rotation: +x,-x-y,+y-z>
>>> print(rot)
+x,-y,-z
>>> print rot.matrix
array([[ 1.,  0.,  0.],
       [ -1.,  1.,  0.],
       [0.,  1.,  -1.]])
exception DoesNotExist
exception MultipleObjectsReturned

Structural

class qmpy.PDF(structure, limit=10)[source]

Container class for a Pair-distribution function.

Attributes:

structure: Structure pairs: dict of (atom1, atom2):[distances] limit: maximum distance

get_pair_distances()[source]

Loops over pairs of atoms that are within radius max_dist of one another. Returns a dict of (atom1, atom2):[list of distances].

class qmpy.XRD(structure=None, measured=False, wavelength=1.5418, min_2th=10, max_2th=60, resolution=0.01)[source]

Container for an X-ray diffraction pattern.

Attributes:
peaks (List):

List of Peak instances.

measured (bool):

True if the XRD is a measured pattern, otherwise False.

min_2th (float):

Minimum 2theta angle allowed. Defaults to 60 degrees.

max_2th (float):

Maximum 2theta angle allowed. Defaults to 10 degrees.

wavelength (float):

X-ray wavelength. Defaults to 1.5418 Ang.

resolution (float):

Minimum 2theta angle the XRD will distinguish between.

get_intensities(bfactors=None, scale=None)[source]

Loops over all peaks calculating intensity.

Keyword Arguments:
bfactors (list)list of B factors for each atomic site. Care must

taken to ensure that the order of B factors agrees with the order of atomic orbits.

scale (float)Scaling factor to multiply the intensities by. If

scale evaluates to False, the intensities will be re-normalized at the end such that the highest peak is 1.

class qmpy.Peak(angle, multiplicity=None, intensity=None, hkl=None, xrd=None, width=None, measured=False)[source]
Attributes:
angle (float):

Peak 2*theta angle in radians.

hkl (list):

HKL indices of the peak.

multiplicity (int):

Number of HKL indices which generate the peak.

lp_factor()[source]

Calculates the Lorentz-polarization factor.

http://reference.iucr.org/dictionary/Lorentz%E2%80%93polarization_correction

thermal_factor(bfactor=1.0)[source]

Calculates the Debye-Waller factor for a peak.

http://en.wikipedia.org/wiki/Debye-Waller_factor