pymatgen.analysis.pourbaix_diagram module

class IonEntry(ion, energy, name=None)[source]

Bases: pymatgen.analysis.phase_diagram.PDEntry

Object similar to PDEntry, but contains an Ion object instead of a Composition object.

Parameters
  • comp – Ion object

  • energy – Energy for composition.

  • name – Optional parameter to name the entry. Defaults to the chemical formula.

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.

as_dict()[source]

Creates a dict of composition, energy, and ion name

classmethod from_dict(d)[source]

Returns an IonEntry object from a dict.

class MultiEntry(entry_list, weights=None)[source]

Bases: pymatgen.analysis.pourbaix_diagram.PourbaixEntry

PourbaixEntry-like object for constructing multi-elemental Pourbaix diagrams.

Initializes a MultiEntry.

Parameters
  • entry_list ([PourbaixEntry]) – List of component PourbaixEntries

  • weights ([float]) – Weights associated with each entry. Default is None

as_dict()[source]

Returns dict which contains Pourbaix Entry data. Note that the pH, voltage, H2O factors are always calculated when constructing a PourbaixEntry object.

classmethod from_dict(d)[source]

Invokes

name

Multientry name, i. e. the name of each entry joined by ‘ + ‘

class PourbaixDiagram(entries, comp_dict=None, conc_dict=None, filter_solids=False, nproc=None)[source]

Bases: monty.json.MSONable

Class to create a Pourbaix diagram from entries

Parameters
  • entries ([PourbaixEntry] or [MultiEntry]) – Entries list containing Solids and Ions or a list of MultiEntries

  • ({str (conc_dict) – float}): Dictionary of compositions, defaults to equal parts of each elements

  • ({str – float}): Dictionary of ion concentrations, defaults to 1e-6 for each element

  • filter_solids (bool) – applying this filter to a pourbaix diagram ensures all included phases are filtered by stability on the compositional phase diagram. This breaks some of the functionality of the analysis, though, so use with caution.

  • nproc (int) – number of processes to generate multientries with in parallel. Defaults to None (serial processing)

all_entries

Return all entries used to generate the pourbaix diagram

as_dict(include_unprocessed_entries=False)[source]

A JSON serializable dict representation of an object.

find_stable_entry(pH, V)[source]

Finds stable entry at a pH,V condition :param pH: pH to find stable entry :type pH: float :param V: V to find stable entry :type V: float

Returns:

classmethod from_dict(d)[source]
get_decomposition_energy(entry, pH, V)[source]

Finds decomposition to most stable entry

Parameters
  • entry (PourbaixEntry) – PourbaixEntry corresponding to compound to find the decomposition for

  • pH (float) – pH at which to find the decomposition

  • V (float) – voltage at which to find the decomposition

Returns

reaction corresponding to the decomposition

get_hull_energy(pH, V)[source]
static get_pourbaix_domains(pourbaix_entries, limits=None)[source]

Returns a set of pourbaix stable domains (i. e. polygons) in pH-V space from a list of pourbaix_entries

This function works by using scipy’s HalfspaceIntersection function to construct all of the 2-D polygons that form the boundaries of the planes corresponding to individual entry gibbs free energies as a function of pH and V. Hyperplanes of the form a*pH + b*V + 1 - g(0, 0) are constructed and supplied to HalfspaceIntersection, which then finds the boundaries of each pourbaix region using the intersection points.

Parameters
  • pourbaix_entries ([PourbaixEntry]) – Pourbaix entries with which to construct stable pourbaix domains

  • limits ([[float]]) – limits in which to do the pourbaix analysis

Returns

[boundary_points]}. The list of boundary points are the sides of the N-1 dim polytope bounding the allowable ph-V range of each entry.

Return type

Returns a dict of the form {entry

static process_multientry(entry_list, prod_comp, coeff_threshold=0.0001)[source]

Static method for finding a multientry based on a list of entries and a product composition. Essentially checks to see if a valid aqueous reaction exists between the entries and the product composition and returns a MultiEntry with weights according to the coefficients if so.

Parameters
  • entry_list ([Entry]) – list of entries from which to create a MultiEntry

  • prod_comp (Composition) – composition constraint for setting weights of MultiEntry

  • coeff_threshold (float) – threshold of stoichiometric coefficients to filter, if weights are lower than this value, the entry is not returned

stable_entries

Returns the stable entries in the Pourbaix diagram.

unprocessed_entries

Return unprocessed entries

unstable_entries

Returns all unstable entries in the Pourbaix diagram

class PourbaixEntry(entry, entry_id=None, concentration=1e-06)[source]

Bases: monty.json.MSONable

An object encompassing all data relevant to a solid or ion in a pourbaix diagram. Each bulk solid/ion has an energy g of the form: e = e0 + 0.0591 log10(conc) - nO mu_H2O + (nH - 2nO) pH + phi (-nH + 2nO + q)

Note that the energies corresponding to the input entries should be formation energies with respect to hydrogen and oxygen gas in order for the pourbaix diagram formalism to work. This may be changed to be more flexible in the future.

Parameters

entry (ComputedEntry/ComputedStructureEntry/PDEntry/IonEntry) – An entry object

as_dict()[source]

Returns dict which contains Pourbaix Entry data. Note that the pH, voltage, H2O factors are always calculated when constructing a PourbaixEntry object.

composition

Returns composition

conc_term

Returns the concentration contribution to the free energy, and should only be present when there are ions in the entry

energy

returns energy

Returns (float): total energy of the pourbaix

entry (at pH, V = 0 vs. SHE)

energy_at_conditions(pH, V)[source]

Get free energy for a given pH and V

Parameters
  • pH (float) – pH at which to evaluate free energy

  • V (float) – voltage at which to evaluate free energy

Returns

free energy at conditions

energy_per_atom

energy per atom of the pourbaix entry

Returns (float): energy per atom

classmethod from_dict(d)[source]

Invokes

nH2O
nPhi
name
normalization_factor

Sum of number of atoms minus the number of H and O in composition

normalized_energy

Returns: energy normalized by number of non H or O atoms, e. g. for Zn2O6, energy / 2 or for AgTe3(OH)3, energy / 4

normalized_energy_at_conditions(pH, V)[source]

Energy at an electrochemical condition, compatible with numpy arrays for pH/V input

Parameters
  • pH (float) – pH at condition

  • V (float) – applied potential at condition

Returns

energy normalized by number of non-O/H atoms at condition

npH
num_atoms

Return number of atoms in current formula. Useful for normalization

class PourbaixPlotter(pourbaix_diagram)[source]

Bases: object

A plotter class for phase diagrams.

Parameters
  • phasediagram – A PhaseDiagram object.

  • show_unstable – Whether unstable phases will be plotted as well as red crosses. Defaults to False.

domain_vertices(entry)[source]

Returns the vertices of the Pourbaix domain.

Parameters

entry – Entry for which domain vertices are desired

Returns

list of vertices

get_pourbaix_plot(limits=None, title='', label_domains=True, plt=None)[source]

Plot Pourbaix diagram.

Parameters
  • limits – 2D list containing limits of the Pourbaix diagram of the form [[xlo, xhi], [ylo, yhi]]

  • title (str) – Title to display on plot

  • label_domains (bool) – whether to label pourbaix domains

  • plt (pyplot) – Pyplot instance for plotting

Returns

plt (pyplot) - matplotlib plot object with pourbaix diagram

plot_entry_stability(entry, pH_range=None, pH_resolution=100, V_range=None, V_resolution=100, e_hull_max=1, cmap='RdYlBu_r', **kwargs)[source]
show(*args, **kwargs)[source]

Shows the pourbaix plot

Parameters
  • *args – args to get_pourbaix_plot

  • **kwargs – kwargs to get_pourbaix_plot

Returns

None

generate_entry_label(entry)[source]

Generates a label for the pourbaix plotter

Parameters

entry (PourbaixEntry or MultiEntry) – entry to get a label for

ion_or_solid_comp_object(formula)[source]

Returns either an ion object or composition object given a formula.

Parameters

formula – String formula. Eg. of ion: NaOH(aq), Na[+]; Eg. of solid: Fe2O3(s), Fe(s), Na2O

Returns

Composition/Ion object

latexify_ion(formula)[source]