pymatgen.analysis.phase_diagram module

This module defines tools to generate and analyze phase diagrams.

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

Bases: 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[source]
as_dict()[source]
Returns:

MSONable dictionary representation of CompoundPhaseDiagram

classmethod from_dict(d)[source]
Parameters:

d (dict) – dictionary representation of CompoundPhaseDiagram

Returns:

CompoundPhaseDiagram

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: 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]
Returns:

MSONable dictionary representation of GrandPotPDEntry

property chemical_energy[source]

The chemical energy term mu*N in the grand potential

Returns:

The chemical energy term mu*N in the grand potential

property composition: Composition[source]

The composition after removing free species

Returns:

Composition

property energy[source]

Returns: The grand potential energy

classmethod from_dict(d)[source]
Parameters:

d (dict) – dictionary representation of GrandPotPDEntry

Returns:

GrandPotPDEntry

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

Bases: 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 entries themselves.

  • computed_data (dict) – A dict containing pre-computed data. This allows PhaseDiagram object to be reconstituted without performing the expensive convex hull computation. The dict is the output from the PhaseDiagram._compute() method and is stored in PhaseDiagram.computed_data when generated for the first time.

as_dict()[source]
Returns:

MSONable dictionary representation of GrandPotentialPhaseDiagram

classmethod from_dict(d)[source]
Parameters:

d (dict) – dictionary representation of GrandPotentialPhaseDiagram

Returns:

GrandPotentialPhaseDiagram

class PDEntry(composition: Composition, energy: float, name: str | None = None, attribute: object = None)[source]

Bases: Entry

An object encompassing all relevant data for phase diagrams.

composition[source]

The composition associated with the PDEntry.

Type:

Composition

energy[source]

The energy associated with the entry.

Type:

float

name[source]

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.

Type:

str

attribute[source]

A arbitrary attribute. Can be used to specify that the entry is a newly found compound, or to specify a particular label for the entry, etc. An attribute can be anything but must be MSONable.

Type:

MSONable

Parameters:
  • composition (Composition) – Composition

  • energy (float) – Energy for composition.

  • name (str) – Optional parameter to name the entry. Defaults to the reduced chemical formula.

  • attribute – Optional attribute of the entry. Must be MSONable.

as_dict()[source]
Returns:

MSONable dictionary representation of PDEntry

property energy: float[source]

Returns: the energy of the entry.

classmethod from_dict(d)[source]
Parameters:

d (dict) – dictionary representation of PDEntry

Returns:

PDEntry

class PDPlotter(phasediagram: PhaseDiagram, show_unstable: float = 0.2, backend: Literal['plotly', 'matplotlib'] = 'plotly', **plotkwargs)[source]

Bases: object

A plotter class for compositional phase diagrams.

Parameters:
  • phasediagram (PhaseDiagram) – PhaseDiagram object.

  • show_unstable (float) – Whether unstable (above the hull) phases will be plotted. If a number > 0 is entered, all phases with e_hull < show_unstable (eV/atom) will be shown.

  • backend ("plotly" | "matplotlib") – Python package used for plotting. Defaults to “plotly”.

  • **plotkwargs (dict) –

    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, plt=None, label_uncertainties=False)[source]
Parameters:
  • label_stable – Whether to label stable compounds.

  • label_unstable – Whether to label unstable compounds.

  • ordering – Ordering of vertices (matplotlib backend only).

  • energy_colormap – Colormap for coloring energy (matplotlib backend only).

  • process_attributes – Whether to process the attributes (matplotlib backend only).

  • plt – Existing plt object if plotting multiple phase diagrams ( matplotlib backend only).

  • label_uncertainties – Whether to add error bars to the hull (plotly backend only). For binaries, this also shades the hull with the uncertainty window.

Returns:

go.Figure (plotly) or matplotlib.pyplot (matplotlib)

property pd_plot_data[source]

Plotting data for phase diagram. Cached for repetitive calls. 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 dict of {coordinatesentry} for each stable node

    in the phase diagram. (Each coordinate can only have one stable phase)

  • unstable_entries is a dict of {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.

plot_element_profile(element, comp, show_label_index=None, xlim=5)[source]

Draw the element profile plot for a composition varying different chemical potential of an element. X value is the negative value of the chemical potential reference to elemental chemical potential. For example, if choose Element(“Li”), X= -(µLi-µLi0), which corresponds to the voltage versus metal anode. Y values represent for the number of element uptake in this composition (unit: per atom). All reactions are printed to help choosing the profile steps you want to show label in the plot.

Parameters:
  • element (Element) – An element of which the chemical potential is considered. It also must be in the phase diagram.

  • comp (Composition) – A composition.

  • show_label_index (list of integers) – The labels for reaction products you want to show in the plot. Default to None (not showing any annotation for reaction products). For the profile steps you want to show the labels, just add it to the show_label_index. The profile step counts from zero. For example, you can set show_label_index=[0, 2, 5] to label profile step 0,2,5.

  • xlim (float) – The max x value. x value is from 0 to xlim. Default to 5 eV.

Returns:

Plot of element profile evolution by varying the chemical potential of an element.

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

Draw the phase diagram using Plotly (or Matplotlib) and show it.

Parameters:
  • *args – Passed to get_plot.

  • **kwargs – Passed to get_plot.

write_image(stream: str | StringIO, image_format: str = 'svg', **kwargs) None[source]

Writes the phase diagram to an image in a stream.

Parameters:
  • stream (str | StringIO) – stream to write to. Can be a file stream or a StringIO stream.

  • image_format (str) – 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 function.

class PatchedPhaseDiagram(entries: Sequence[PDEntry] | set[PDEntry], elements: Sequence[Element] | None = None, keep_all_spaces: bool = False, verbose: bool = False)[source]

Bases: PhaseDiagram

Computing the Convex Hull of a large set of data in multiple dimensions is highly expensive. This class acts to breakdown large chemical spaces into smaller chemical spaces which can be computed much more quickly due to having both reduced dimensionality and data set sizes.

subspaces ({str

{Element, }}): Dictionary of the sets of elements for each of the PhaseDiagrams within the PatchedPhaseDiagram.

pds ({str

PhaseDiagram}): Dictionary of PhaseDiagrams within the PatchedPhaseDiagram.

all_entries[source]

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.

Type:

list[PDEntry]

min_entries[source]

List of the lowest energy entries for each composition in the data provided for Phase Diagram construction.

Type:

list[PDEntry]

el_refs[source]

List of elemental references for the phase diagrams. These are entries corresponding to the lowest energy element entries for simple compositional phase diagrams.

Type:

list[PDEntry]

elements[source]

List of elements in the phase diagram.

Type:

list[Element]

Parameters:
  • entries (list[PDEntry]) – A list of PDEntry-like objects having an energy, energy_per_atom and composition.

  • elements (list[Element], optional) – Optional list of elements in the phase diagram. If set to None, the elements are determined from the entries themselves and are sorted alphabetically. If specified, element ordering (e.g. for pd coordinates) is preserved.

  • keep_all_spaces (bool) – Boolean control on whether to keep chemical spaces that are subspaces of other spaces.

  • verbose (bool) – Whether to show progress bar during convex hull construction.

as_dict() dict[str, Any][source]
Returns:

MSONable dictionary representation of PatchedPhaseDiagram

Return type:

dict[str, Any]

classmethod from_dict(d)[source]
Parameters:

d (dict) – dictionary representation of PatchedPhaseDiagram

Returns:

PatchedPhaseDiagram

get_all_chempots()[source]

Not Implemented - See PhaseDiagram

get_chempot_range_map()[source]

Not Implemented - See PhaseDiagram

get_chempot_range_stability_phase()[source]

Not Implemented - See PhaseDiagram

get_composition_chempots()[source]

Not Implemented - See PhaseDiagram

get_critical_compositions()[source]

Not Implemented - See PhaseDiagram

get_decomp_and_e_above_hull(entry: PDEntry, allow_negative: bool = False, check_stable: bool = False, on_error: Literal['raise', 'warn', 'ignore'] = 'raise') tuple[dict[PDEntry, float], float] | tuple[None, None][source]

Same as method on parent class PhaseDiagram except check_stable defaults to False for speed. See https://github.com/materialsproject/pymatgen/issues/2840 for details.

get_decomposition(comp: Composition) dict[pymatgen.analysis.phase_diagram.PDEntry, float][source]

See PhaseDiagram

Parameters:

comp (Composition) – A composition

Returns:

amount} where amount is the amount of the fractional composition.

Return type:

Decomposition as a dict of {PDEntry

get_element_profile()[source]

Not Implemented - See PhaseDiagram

get_equilibrium_reaction_energy(entry: Entry) float[source]

See PhaseDiagram

NOTE this is only approximately the same as the what we would get from PhaseDiagram as we make use of the slsqp approach inside get_phase_separation_energy().

Parameters:

entry (PDEntry) – A PDEntry like object

Returns:

Equilibrium reaction energy of entry. Stable entries should have equilibrium reaction energy <= 0. The energy is given per atom.

get_pd_for_entry(entry: Entry | Composition) PhaseDiagram[source]

Get the possible phase diagrams for an entry

Parameters:

entry (PDEntry | Composition) – A PDEntry or Composition-like object

Returns:

phase diagram that the entry is part of

Return type:

PhaseDiagram

get_transition_chempots()[source]

Not Implemented - See PhaseDiagram

getmu_vertices_stability_phase()[source]

Not Implemented - See PhaseDiagram

class PhaseDiagram(entries: Sequence[PDEntry] | set[PDEntry], elements: Sequence[Element] = (), *, computed_data: dict[str, Any] | None = None)[source]

Bases: MSONable

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

      1. 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

      1. 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

dim[source]

The dimensionality of the phase diagram.

Type:

int

elements[source]

Elements in the phase diagram.

el_refs[source]

List of elemental references for the phase diagrams. These are entries corresponding to the lowest energy element entries for simple compositional phase diagrams.

all_entries[source]

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.

qhull_entries[source]

Actual entries used in convex hull. Excludes all positive formation energy entries.

qhull_data[source]

Data used in the convex hull operation. This is essentially a matrix of composition data and energy per atom values created from qhull_entries.

facets[source]

Facets of the phase diagram in the form of [[1,2,3],[4,5,6]…]. For a ternary, it is the indices (references to qhull_entries and qhull_data) for the vertices of the phase triangles. Similarly extended to higher D simplices for higher dimensions.

simplices[source]

The simplices of the phase diagram as a list of np.ndarray, i.e., the list of stable compositional coordinates in the phase diagram.

Parameters:
  • entries (list[PDEntry]) – A list of PDEntry-like objects having an energy, energy_per_atom and composition.

  • elements (list[Element]) – Optional list of elements in the phase diagram. If set to None, the elements are determined from the entries themselves and are sorted alphabetically. If specified, element ordering (e.g. for pd coordinates) is preserved.

  • computed_data (dict) – A dict containing pre-computed data. This allows PhaseDiagram object to be reconstituted without performing the expensive convex hull computation. The dict is the output from the PhaseDiagram._compute() method and is stored in PhaseDiagram.computed_data when generated for the first time.

property all_entries_hulldata[source]

Returns: The actual ndarray used to construct the convex hull.

as_dict()[source]
Returns:

MSONable dictionary representation of PhaseDiagram

formation_energy_tol = 1e-11[source]
classmethod from_dict(d: dict[str, Any]) PhaseDiagram[source]
Parameters:

d (dict) – dictionary representation of PhaseDiagram

Returns:

PhaseDiagram

get_all_chempots(comp)[source]

Get chemical potentials at a given composition.

Parameters:

comp (Composition) – Composition

Returns:

Chemical potentials.

get_chempot_range_map(elements: Sequence[Element], referenced: bool = True, joggle: bool = True) dict[pymatgen.core.periodic_table.Element, list[pymatgen.util.coord.Simplex]][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 (bool) – 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 corresponding 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 the chemical potentials for all elements at a given composition.

Parameters:

comp (Composition) – Composition

Returns:

Dictionary of chemical potentials.

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.

Parameters:
  • comp1 (Composition) – First composition to define the tieline

  • comp2 (Composition) – Second composition to define the tieline

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: PDEntry, allow_negative: bool = False, check_stable: bool = True, on_error: Literal['raise', 'warn', 'ignore'] = 'raise') tuple[dict[PDEntry, float], float] | tuple[None, None][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 (PDEntry) – A PDEntry like object

  • allow_negative (bool) – Whether to allow negative e_above_hulls. Used to calculate equilibrium reaction energies. Defaults to False.

  • check_stable (bool) – Whether to first check whether an entry is stable. In normal circumstances, this is the faster option since checking for stable entries is relatively fast. However, if you have a huge proportion of unstable entries, then this check can slow things down. You should then set this to False.

  • on_error ('raise' | 'warn' | 'ignore') – What to do if no valid decomposition was found. ‘raise’ will throw ValueError. ‘warn’ will print return (None, None). ‘ignore’ just returns (None, None). Defaults to ‘raise’.

Raises:

ValueError – If no valid decomposition exists in this phase diagram for given entry.

Returns:

(decomp, energy_above_hull). The decomposition is provided

as a dict of {PDEntry: amount} where amount is the amount of the fractional composition. Stable entries should have energy above convex hull of 0. The energy is given per atom.

get_decomp_and_hull_energy_per_atom(comp: Composition) tuple[dict[pymatgen.analysis.phase_diagram.PDEntry, float], float][source]
Parameters:

comp (Composition) – Input composition

Returns:

Energy of lowest energy equilibrium at desired composition per atom

get_decomp_and_phase_separation_energy(entry: PDEntry, space_limit: int = 200, stable_only: bool = False, tols: Sequence[float] = (1e-08,), maxiter: int = 1000, **kwargs: Any) tuple[dict[PDEntry, float], float] | tuple[None, None][source]

Provides the combination of entries in the PhaseDiagram that gives the lowest formation enthalpy with the same composition as the given entry excluding entries with the same composition and the energy difference per atom between the given entry and the energy of the combination found.

For unstable entries that are not polymorphs of stable entries (or completely novel entries) this is simply the energy above (or below) the convex hull.

For entries with the same composition as one of the stable entries in the phase diagram setting stable_only to False (Default) allows for entries not previously on the convex hull to be considered in the combination. In this case the energy returned is what is referred to as the decomposition enthalpy in:

  1. Bartel, C., Trewartha, A., Wang, Q., Dunn, A., Jain, A., Ceder, G.,

    A critical examination of compound stability predictions from machine-learned formation energies, npj Computational Materials 6, 97 (2020)

For stable entries setting stable_only to True returns the same energy as get_equilibrium_reaction_energy. This function is based on a constrained optimization rather than recalculation of the convex hull making it algorithmically cheaper. However, if tol is too loose there is potential for this algorithm to converge to a different solution.

Parameters:
  • entry (PDEntry) – A PDEntry like object.

  • space_limit (int) – The maximum number of competing entries to consider before calculating a second convex hull to reducing the complexity of the optimization.

  • stable_only (bool) – Only use stable materials as competing entries.

  • tols (list[float]) – Tolerances for convergence of the SLSQP optimization when finding the equilibrium reaction. Tighter tolerances tested first.

  • maxiter (int) – The maximum number of iterations of the SLSQP optimizer when finding the equilibrium reaction.

  • **kwargs – Passed to get_decomp_and_e_above_hull.

Returns:

(decomp, energy). The decomposition is given as a dict of {PDEntry, amount} for all entries in the decomp reaction where amount is the amount of the fractional composition. The phase separation energy is given per atom.

get_decomposition(comp: Composition) dict[pymatgen.analysis.phase_diagram.PDEntry, float][source]

Provides the decomposition at a particular composition.

Parameters:

comp (Composition) – A composition

Returns:

amount} where amount is the amount of the fractional composition.

Return type:

Decomposition as a dict of {PDEntry

get_e_above_hull(entry: PDEntry, **kwargs: Any) float | None[source]

Provides the energy above convex hull for an entry

Parameters:
  • entry (PDEntry) – A PDEntry like object.

  • **kwargs – Passed to get_decomp_and_e_above_hull().

Returns:

Energy above convex hull of entry. Stable entries should have

energy above hull of 0. The energy is given per atom.

Return type:

float | None

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: PDEntry) float | None[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 (PDEntry) – A PDEntry like object

Returns:

Equilibrium reaction energy of entry. Stable entries should have

equilibrium reaction energy <= 0. The energy is given per atom.

Return type:

float | None

get_form_energy(entry: PDEntry) float[source]

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

Parameters:

entry (PDEntry) – A PDEntry-like object.

Returns:

Formation energy from the elemental references.

Return type:

float

get_form_energy_per_atom(entry: PDEntry) float[source]

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

Parameters:

entry (PDEntry) – An PDEntry-like object

Returns:

Formation energy per atom from the elemental references.

get_hull_energy(comp: Composition) float[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_hull_energy_per_atom(comp: Composition, **kwargs) float[source]
Parameters:

comp (Composition) – Input composition

Returns:

Energy of lowest energy equilibrium at desired composition.

get_phase_separation_energy(entry, **kwargs)[source]

Provides the energy to the convex hull for the given entry. For stable entries already in the phase diagram the algorithm provides the phase separation energy which is referred to as the decomposition enthalpy in:

  1. Bartel, C., Trewartha, A., Wang, Q., Dunn, A., Jain, A., Ceder, G.,

    A critical examination of compound stability predictions from machine-learned formation energies, npj Computational Materials 6, 97 (2020)

Parameters:
  • entry (PDEntry) – A PDEntry like object

  • **kwargs

    Keyword args passed to get_decomp_and_decomp_energy space_limit (int): The maximum number of competing entries to consider. stable_only (bool): Only use stable materials as competing entries tol (float): The tolerance for convergence of the SLSQP optimization

    when finding the equilibrium reaction.

    maxiter (int): The maximum number of iterations of the SLSQP optimizer

    when finding the equilibrium reaction.

Returns:

phase separation energy per atom of entry. Stable entries should have energies <= 0, Stable elemental entries should have energies = 0 and unstable entries should have energies > 0. Entries that have the same composition as a stable energy may have positive or negative phase separation energies depending on their own energy.

get_reference_energy_per_atom(comp: Composition) float[source]
Parameters:

comp (Composition) – Input composition

Returns:

Reference energy of the terminal species at a given composition.

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 the stable phase at the target composition

  • 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[source]
pd_coords(comp: Composition) ndarray[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.

Parameters:

comp (Composition) – A composition

Returns:

The coordinates for a given composition in the PhaseDiagram’s basis

property stable_entries: set[pymatgen.entries.Entry][source]

Returns: set[Entry]: of stable entries in the phase diagram.

property unstable_entries: set[pymatgen.entries.Entry][source]

Returns: set[Entry]: unstable entries in the phase diagram. Includes positive formation energy entries.

exception PhaseDiagramError[source]

Bases: Exception

An exception class for Phase Diagram generation.

class ReactionDiagram(entry1, entry2, all_entries, tol: float = 0.0001, float_fmt='%.4f')[source]

Bases: object

Analyzes the possible reactions between a pair of compounds, e.g., an electrolyte and an electrode.

Parameters:
  • entry1 (ComputedEntry) – Entry for 1st component. Note that corrections, if any, must already be pre-applied. This is to give flexibility for different kinds of corrections, e.g., if a particular entry is fitted to an experimental data (such as EC molecule).

  • entry2 (ComputedEntry) – Entry for 2nd component. Note that corrections must already be pre-applied. This is to give flexibility for different kinds of corrections, e.g., if a particular entry is fitted to an experimental data (such as EC molecule).

  • all_entries ([ComputedEntry]) – All other entries to be considered in the analysis. Note that corrections, if any, must already be pre-applied.

  • tol (float) – Tolerance to be used to determine validity of reaction.

  • float_fmt (str) – Formatting string to be applied to all floats. Determines number of decimal places in reaction string.

get_compound_pd()[source]

Get the CompoundPhaseDiagram object, which can then be used for plotting.

Returns:

CompoundPhaseDiagram

class TransformedPDEntry(entry, sp_mapping, name=None)[source]

Bases: PDEntry

This class represents 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:
  • entry (PDEntry) – Original entry to be transformed.

  • ({Composition (sp_mapping) – DummySpecies}): dictionary mapping Terminal Compositions to Dummy Species

amount_tol = 1e-05[source]
as_dict()[source]
Returns:

MSONable dictionary representation of TransformedPDEntry

property composition: Composition[source]

The composition in the dummy species space

Returns:

Composition

classmethod from_dict(d)[source]
Parameters:

d (dict) – dictionary representation of TransformedPDEntry

Returns:

TransformedPDEntry

exception TransformedPDEntryError[source]

Bases: Exception

An exception class for TransformedPDEntry.

get_facets(qhull_data: ArrayLike, joggle: bool = False) ConvexHull[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 (bool) – 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:

coord – 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:

coord – 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