pymatgen.analysis.phase_diagram module

class CompoundPhaseDiagram(entries, terminal_compositions, normalize_terminal_compositions=True)[source]

Bases: pymatgen.analysis.phase_diagram.PhaseDiagram

Generates phase diagrams from compounds as terminations instead of elements.

Initializes a CompoundPhaseDiagram.

Parameters:
  • entries ([PDEntry]) – Sequence of input entries. For example, if you want a Li2O-P2O5 phase diagram, you might have all Li-P-O entries as an input.
  • terminal_compositions ([Composition]) – Terminal compositions of phase space. In the Li2O-P2O5 example, these will be the Li2O and P2O5 compositions.
  • normalize_terminal_compositions (bool) – Whether to normalize the terminal compositions to a per atom basis. If normalized, the energy above hulls will be consistent for comparison across systems. Non-normalized terminals are more intuitive in terms of compositional breakdowns.
amount_tol = 1e-05
as_dict()[source]
classmethod from_dict(d)[source]
transform_entries(entries, terminal_compositions)[source]

Method to transform all entries to the composition coordinate in the terminal compositions. If the entry does not fall within the space defined by the terminal compositions, they are excluded. For example, Li3PO4 is mapped into a Li2O:1.5, P2O5:0.5 composition. The terminal compositions are represented by DummySpecies.

Parameters:
  • entries – Sequence of all input entries
  • terminal_compositions – Terminal compositions of phase space.
Returns:

Sequence of TransformedPDEntries falling within the phase space.

class GrandPotPDEntry(entry, chempots, name=None)[source]

Bases: pymatgen.analysis.phase_diagram.PDEntry

A grand potential pd entry object encompassing all relevant data for phase diagrams. Chemical potentials are given as a element-chemical potential dict.

Parameters:
  • entry – A PDEntry-like object.
  • chempots – Chemical potential specification as {Element: float}.
  • name – Optional parameter to name the entry. Defaults to the reduced chemical formula of the original entry.
as_dict()[source]
classmethod from_dict(d)[source]
is_element

True if the entry is an element.

class GrandPotentialPhaseDiagram(entries, chempots, elements=None)[source]

Bases: pymatgen.analysis.phase_diagram.PhaseDiagram

A class representing a Grand potential phase diagram. Grand potential phase diagrams are essentially phase diagrams that are open to one or more components. To construct such phase diagrams, the relevant free energy is the grand potential, which can be written as the Legendre transform of the Gibbs free energy as follows

Grand potential = G - u_X N_X

The algorithm is based on the work in the following papers:

  1. S. P. Ong, L. Wang, B. Kang, and G. Ceder, Li-Fe-P-O2 Phase Diagram from First Principles Calculations. Chem. Mater., 2008, 20(5), 1798-1807. doi:10.1021/cm702327g
  2. S. P. Ong, A. Jain, G. Hautier, B. Kang, G. Ceder, Thermal stabilities of delithiated olivine MPO4 (M=Fe, Mn) cathodes investigated using first principles calculations. Electrochem. Comm., 2010, 12(3), 427-430. doi:10.1016/j.elecom.2010.01.010

Standard constructor for grand potential phase diagram.

Parameters:
  • entries ([PDEntry]) – A list of PDEntry-like objects having an energy, energy_per_atom and composition.
  • {Element (chempots) – float}: Specify the chemical potentials of the open elements.
  • elements ([Element]) – Optional list of elements in the phase diagram. If set to None, the elements are determined from the the entries themselves.
as_dict()[source]
classmethod from_dict(d)[source]
class PDAnalyzer(*args, **kwargs)[source]

Bases: object

class PDEntry(composition, energy, name=None, attribute=None)[source]

Bases: monty.json.MSONable

An object encompassing all relevant data for phase diagrams.

name

A name for the entry. This is the string shown in the phase diagrams. By default, this is the reduced formula for the composition, but can be set to some other string for display purposes.

Parameters:
  • comp – Composition as a pymatgen.core.structure.Composition
  • energy – Energy for composition.
  • name – Optional parameter to name the entry. Defaults to the reduced chemical formula.
  • attribute – Optional attribute of the entry. This can be used to specify that the entry is a newly found compound, or to specify a particular label for the entry, or else … Used for further analysis and plotting purposes. An attribute can be anything but must be MSONable.
as_dict()[source]
energy_per_atom

Returns the final energy per atom.

static from_csv(filename)[source]

Imports PDEntries from a csv.

Parameters:filename – Filename to import from.
Returns:List of Elements, List of PDEntries
classmethod from_dict(d)[source]
is_element

True if the entry is an element.

static to_csv(filename, entries, latexify_names=False)[source]

Exports PDEntries to a csv

Parameters:
  • filename – Filename to write to.
  • entries – PDEntries to export.
  • latexify_names – Format entry names to be LaTex compatible, e.g., Li_{2}O
class PDEntryIO[source]

Bases: object

static from_csv(*args, **kwargs)[source]
static to_csv(*args, **kwargs)[source]
class PDPlotter(phasediagram, show_unstable=0, **plotkwargs)[source]

Bases: object

A plotter class for phase diagrams.

Parameters:
  • phasediagram – PhaseDiagram object.
  • show_unstable (float) – Whether unstable phases will be plotted as well as red crosses. If a number > 0 is entered, all phases with ehull < show_unstable will be shown.
  • **plotkwargs

    Keyword args passed to matplotlib.pyplot.plot. Can be used to customize markers etc. If not set, the default is {

    ”markerfacecolor”: (0.2157, 0.4941, 0.7216), “markersize”: 10, “linewidth”: 3

    }

get_chempot_range_map_plot(elements, referenced=True)[source]

Returns a plot of the chemical potential range _map. Currently works only for 3-component PDs.

Parameters:
  • elements – Sequence of elements to be considered as independent variables. E.g., if you want to show the stability ranges of all Li-Co-O phases wrt to uLi and uO, you will supply [Element(“Li”), Element(“O”)]
  • referenced – if True, gives the results with a reference being the energy of the elemental phase. If False, gives absolute values.
Returns:

A matplotlib plot object.

get_contour_pd_plot()[source]

Plot a contour phase diagram plot, where phase triangles are colored according to degree of instability by interpolation. Currently only works for 3-component phase diagrams.

Returns:A matplotlib plot object.
get_plot(label_stable=True, label_unstable=True, ordering=None, energy_colormap=None, process_attributes=False)[source]
pd_plot_data

Plot data for phase diagram. 2-comp - Full hull with energies 3/4-comp - Projection into 2D or 3D Gibbs triangle.

Returns:
  • lines is a list of list of coordinates for lines in the PD.
  • stable_entries is a {coordinate : entry} for each stable node

in the phase diagram. (Each coordinate can only have one stable phase) - unstable_entries is a {entry: coordinates} for all unstable nodes in the phase diagram.

Return type:(lines, stable_entries, unstable_entries)
plot_chempot_range_map(elements, referenced=True)[source]

Plot the chemical potential range _map. Currently works only for 3-component PDs.

Parameters:
  • elements – Sequence of elements to be considered as independent variables. E.g., if you want to show the stability ranges of all Li-Co-O phases wrt to uLi and uO, you will supply [Element(“Li”), Element(“O”)]
  • referenced – if True, gives the results with a reference being the energy of the elemental phase. If False, gives absolute values.
show(*args, **kwargs)[source]

Draws the phase diagram using Matplotlib and show it.

Parameters:
  • *args – Passed to get_plot.
  • **kwargs – Passed to get_plot.
write_image(stream, image_format='svg', **kwargs)[source]

Writes the phase diagram to an image in a stream.

Parameters:
  • stream – stream to write to. Can be a file stream or a StringIO stream.
  • image_format – format for image. Can be any of matplotlib supported formats. Defaults to svg for best results for vector graphics.
  • **kwargs – Pass through to get_plot functino.
class PhaseDiagram(entries, elements=None)[source]

Bases: monty.json.MSONable

Simple phase diagram class taking in elements and entries as inputs. The algorithm is based on the work in the following papers:

  1. S. P. Ong, L. Wang, B. Kang, and G. Ceder, Li-Fe-P-O2 Phase Diagram from First Principles Calculations. Chem. Mater., 2008, 20(5), 1798-1807. doi:10.1021/cm702327g
  2. S. P. Ong, A. Jain, G. Hautier, B. Kang, G. Ceder, Thermal stabilities of delithiated olivine MPO4 (M=Fe, Mn) cathodes investigated using first principles calculations. Electrochem. Comm., 2010, 12(3), 427-430. doi:10.1016/j.elecom.2010.01.010

..attribute: all_entries

All entries provided for Phase Diagram construction. Note that this does not mean that all these entries are actually used in the phase diagram. For example, this includes the positive formation energy entries that are filtered out before Phase Diagram construction.

Standard constructor for phase diagram.

Parameters:
  • entries ([PDEntry]) – A list of PDEntry-like objects having an energy, energy_per_atom and composition.
  • elements ([Element]) – Optional list of elements in the phase diagram. If set to None, the elements are determined from the the entries themselves.
all_entries_hulldata
as_dict()[source]
formation_energy_tol = 1e-11
classmethod from_dict(d)[source]
get_chempot_range_map(elements, referenced=True, joggle=True)[source]

Returns a chemical potential range map for each stable entry.

Parameters:
  • elements – Sequence of elements to be considered as independent variables. E.g., if you want to show the stability ranges of all Li-Co-O phases wrt to uLi and uO, you will supply [Element(“Li”), Element(“O”)]
  • referenced – If True, gives the results with a reference being the energy of the elemental phase. If False, gives absolute values.
  • joggle (boolean) – Whether to joggle the input to avoid precision errors.
Returns:

[simplices]}. The list of simplices are the sides of the N-1 dim polytope bounding the allowable chemical potential range of each entry.

Return type:

Returns a dict of the form {entry

get_chempot_range_stability_phase(target_comp, open_elt)[source]

returns a set of chemical potentials correspoding to the max and min chemical potential of the open element for a given composition. It is quite common to have for instance a ternary oxide (e.g., ABO3) for which you want to know what are the A and B chemical potential leading to the highest and lowest oxygen chemical potential (reducing and oxidizing conditions). This is useful for defect computations.

Parameters:
  • target_comp – A Composition object
  • open_elt – Element that you want to constrain to be max or min
Returns:

(mu_min,mu_max)}: Chemical potentials are given in “absolute” values (i.e., not referenced to 0)

Return type:

{Element

get_composition_chempots(comp)[source]
get_critical_compositions(comp1, comp2)[source]

Get the critical compositions along the tieline between two compositions. I.e. where the decomposition products change. The endpoints are also returned. :param comp1, comp2: compositions that define the tieline :type comp1, comp2: Composition

Returns:
list of critical compositions. All are of
the form x * comp1 + (1-x) * comp2
Return type:[(Composition)]
get_decomp_and_e_above_hull(entry, allow_negative=False)[source]

Provides the decomposition and energy above convex hull for an entry. Due to caching, can be much faster if entries with the same composition are processed together.

Parameters:
  • entry – A PDEntry like object
  • allow_negative – Whether to allow negative e_above_hulls. Used to calculate equilibrium reaction energies. Defaults to False.
Returns:

(decomp, energy above convex hull) Stable entries should have energy above hull of 0. The decomposition is provided as a dict of {Entry: amount}.

get_decomposition(comp)[source]

Provides the decomposition at a particular composition.

Parameters:comp – A composition
Returns:amount}
Return type:Decomposition as a dict of {Entry
get_e_above_hull(entry)[source]

Provides the energy above convex hull for an entry

Parameters:entry – A PDEntry like object
Returns:Energy above convex hull of entry. Stable entries should have energy above hull of 0.
get_element_profile(element, comp, comp_tol=1e-05)[source]

Provides the element evolution data for a composition. For example, can be used to analyze Li conversion voltages by varying uLi and looking at the phases formed. Also can be used to analyze O2 evolution by varying uO2.

Parameters:
  • element – An element. Must be in the phase diagram.
  • comp – A Composition
  • comp_tol – The tolerance to use when calculating decompositions. Phases with amounts less than this tolerance are excluded. Defaults to 1e-5.
Returns:

[ {‘chempot’: -10.487582010000001, ‘evolution’: -2.0, ‘reaction’: Reaction Object], …]

Return type:

Evolution data as a list of dictionaries of the following format

get_equilibrium_reaction_energy(entry)[source]

Provides the reaction energy of a stable entry from the neighboring equilibrium stable entries (also known as the inverse distance to hull).

Parameters:entry – A PDEntry like object
Returns:Equilibrium reaction energy of entry. Stable entries should have equilibrium reaction energy <= 0.
get_facet_chempots(*args, **kwargs)
get_form_energy(entry)[source]

Returns the formation energy for an entry (NOT normalized) from the elemental references.

Parameters:entry – A PDEntry-like object.
Returns:Formation energy from the elemental references.
get_form_energy_per_atom(entry)[source]

Returns the formation energy per atom for an entry from the elemental references.

Parameters:entry – An PDEntry-like object
Returns:Formation energy per atom from the elemental references.
get_hull_energy(comp)[source]
Parameters:comp (Composition) – Input composition
Returns:Energy of lowest energy equilibrium at desired composition. Not normalized by atoms, i.e. E(Li4O2) = 2 * E(Li2O)
get_transition_chempots(element)[source]

Get the critical chemical potentials for an element in the Phase Diagram.

Parameters:element – An element. Has to be in the PD in the first place.
Returns:A sorted sequence of critical chemical potentials, from less negative to more negative.
getmu_vertices_stability_phase(target_comp, dep_elt, tol_en=0.01)[source]

returns a set of chemical potentials corresponding to the vertices of the simplex in the chemical potential phase diagram. The simplex is built using all elements in the target_composition except dep_elt. The chemical potential of dep_elt is computed from the target composition energy. This method is useful to get the limiting conditions for defects computations for instance.

Parameters:
  • target_comp – A Composition object
  • dep_elt – the element for which the chemical potential is computed from the energy of
  • stable phase at the target composition (the) –
  • tol_en – a tolerance on the energy to set
Returns:

mu}]: An array of conditions on simplex vertices for which each element has a chemical potential set to a given value. “absolute” values (i.e., not referenced to element energies)

Return type:

[{Element

numerical_tol = 1e-08
pd_coords(comp)[source]

The phase diagram is generated in a reduced dimensional space (n_elements - 1). This function returns the coordinates in that space. These coordinates are compatible with the stored simplex objects.

stable_entries

Returns the stable entries in the phase diagram.

unstable_entries

Entries that are unstable in the phase diagram. Includes positive formation energy entries.

exception PhaseDiagramError[source]

Bases: Exception

An exception class for Phase Diagram generation.

class TransformedPDEntry(comp, original_entry)[source]

Bases: pymatgen.analysis.phase_diagram.PDEntry

This class repesents a TransformedPDEntry, which allows for a PDEntry to be transformed to a different composition coordinate space. It is used in the construction of phase diagrams that do not have elements as the terminal compositions.

Parameters:
  • comp – Transformed composition as a Composition.
  • energy – Energy for composition.
  • original_entry – Original entry that this entry arose from.
as_dict()[source]
classmethod from_dict(d)[source]
get_facets(qhull_data, joggle=False)[source]

Get the simplex facets for the Convex hull.

Parameters:
  • qhull_data (np.ndarray) – The data from which to construct the convex hull as a Nxd array (N being number of data points and d being the dimension)
  • joggle (boolean) – Whether to joggle the input to avoid precision errors.
Returns:

List of simplices of the Convex Hull.

order_phase_diagram(lines, stable_entries, unstable_entries, ordering)[source]

Orders the entries (their coordinates) in a phase diagram plot according to the user specified ordering. Ordering should be given as [‘Up’, ‘Left’, ‘Right’], where Up, Left and Right are the names of the entries in the upper, left and right corners of the triangle respectively.

Parameters:
  • lines – list of list of coordinates for lines in the PD.
  • stable_entries – {coordinate : entry} for each stable node in the phase diagram. (Each coordinate can only have one stable phase)
  • unstable_entries – {entry: coordinates} for all unstable nodes in the phase diagram.
  • ordering – Ordering of the phase diagram, given as a list [‘Up’, ‘Left’,’Right’]
Returns:

  • newlines is a list of list of coordinates for lines in the PD.
  • newstable_entries is a {coordinate : entry} for each stable node

in the phase diagram. (Each coordinate can only have one stable phase) - newunstable_entries is a {entry: coordinates} for all unstable nodes in the phase diagram.

Return type:

(newlines, newstable_entries, newunstable_entries)

tet_coord(coord)[source]

Convert a 3D coordinate into a tetrahedron based coordinate system for a prettier phase diagram.

Parameters:coordinate – coordinate used in the convex hull computation.
Returns:coordinates in a tetrahedron-based coordinate system.
triangular_coord(coord)[source]

Convert a 2D coordinate into a triangle-based coordinate system for a prettier phase diagram.

Parameters:coordinate – coordinate used in the convex hull computation.
Returns:coordinates in a triangular-based coordinate system.
uniquelines(q)[source]

Given all the facets, convert it into a set of unique lines. Specifically used for converting convex hull facets into line pairs of coordinates.

Parameters:q – A 2-dim sequence, where each row represents a facet. E.g., [[1,2,3],[3,6,7],…]
Returns:A set of tuple of lines. E.g., ((1,2), (1,3), (2,3), ….)
Return type:setoflines