pymatgen.phasediagram.maker module

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

Bases: pymatgen.phasediagram.maker.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 GrandPotentialPhaseDiagram(entries, chempots, elements=None)[source]

Bases: pymatgen.phasediagram.maker.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 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_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.
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.

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.